xref: /petsc/src/snes/interface/snes.c (revision 29bbc08cd5461c366aba6645924e8ff42acd1de0)
1*29bbc08cSBarry Smith /*$Id: snes.c,v 1.221 2000/09/22 20:45:56 bsmith Exp bsmith $*/
29b94acceSBarry Smith 
3e090d566SSatish Balay #include "src/snes/snesimpl.h"      /*I "petscsnes.h"  I*/
49b94acceSBarry Smith 
54c49b128SBarry Smith PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
6488ecbafSBarry Smith FList      SNESList = 0;
782bf6240SBarry Smith 
85615d1e5SSatish Balay #undef __FUNC__
9b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESView"
107e2c5f70SBarry Smith /*@C
119b94acceSBarry Smith    SNESView - Prints the SNES data structure.
129b94acceSBarry Smith 
134c49b128SBarry Smith    Collective on SNES
14fee21e36SBarry Smith 
15c7afd0dbSLois Curfman McInnes    Input Parameters:
16c7afd0dbSLois Curfman McInnes +  SNES - the SNES context
17c7afd0dbSLois Curfman McInnes -  viewer - visualization context
18c7afd0dbSLois Curfman McInnes 
199b94acceSBarry Smith    Options Database Key:
20c8a8ba5cSLois Curfman McInnes .  -snes_view - Calls SNESView() at end of SNESSolve()
219b94acceSBarry Smith 
229b94acceSBarry Smith    Notes:
239b94acceSBarry Smith    The available visualization contexts include
24c7afd0dbSLois Curfman McInnes +     VIEWER_STDOUT_SELF - standard output (default)
25c7afd0dbSLois Curfman McInnes -     VIEWER_STDOUT_WORLD - synchronized standard
26c8a8ba5cSLois Curfman McInnes          output where only the first processor opens
27c8a8ba5cSLois Curfman McInnes          the file.  All other processors send their
28c8a8ba5cSLois Curfman McInnes          data to the first processor to print.
299b94acceSBarry Smith 
303e081fefSLois Curfman McInnes    The user can open an alternative visualization context with
3177ed5343SBarry Smith    ViewerASCIIOpen() - output to a specified file.
329b94acceSBarry Smith 
3336851e7fSLois Curfman McInnes    Level: beginner
3436851e7fSLois Curfman McInnes 
359b94acceSBarry Smith .keywords: SNES, view
369b94acceSBarry Smith 
3777ed5343SBarry Smith .seealso: ViewerASCIIOpen()
389b94acceSBarry Smith @*/
399b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer)
409b94acceSBarry Smith {
419b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
429b94acceSBarry Smith   int                 ierr;
439b94acceSBarry Smith   SLES                sles;
44454a90a3SBarry Smith   char                *type;
456831982aSBarry Smith   PetscTruth          isascii,isstring;
469b94acceSBarry Smith 
473a40ed3dSBarry Smith   PetscFunctionBegin;
4874679c65SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
493eda8832SBarry Smith   if (!viewer) viewer = VIEWER_STDOUT_(snes->comm);
500f5bd95cSBarry Smith   PetscValidHeaderSpecific(viewer,VIEWER_COOKIE);
516831982aSBarry Smith   PetscCheckSameComm(snes,viewer);
5274679c65SBarry Smith 
536831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr);
546831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,STRING_VIEWER,&isstring);CHKERRQ(ierr);
550f5bd95cSBarry Smith   if (isascii) {
56d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
57454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
58454a90a3SBarry Smith     if (type) {
59454a90a3SBarry Smith       ierr = ViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
60184914b5SBarry Smith     } else {
61454a90a3SBarry Smith       ierr = ViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
62184914b5SBarry Smith     }
630ef38995SBarry Smith     if (snes->view) {
640ef38995SBarry Smith       ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
650ef38995SBarry Smith       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
660ef38995SBarry Smith       ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
670ef38995SBarry Smith     }
68d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
69d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n",
70d132466eSBarry Smith                  snes->rtol,snes->atol,snes->trunctol,snes->xtol);CHKERRQ(ierr);
71d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr);
72d132466eSBarry Smith     ierr = ViewerASCIIPrintf(viewer,"  total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr);
730ef38995SBarry Smith     if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
74d132466eSBarry Smith       ierr = ViewerASCIIPrintf(viewer,"  min function tolerance=%g\n",snes->fmin);CHKERRQ(ierr);
750ef38995SBarry Smith     }
769b94acceSBarry Smith     if (snes->ksp_ewconv) {
779b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
789b94acceSBarry Smith       if (kctx) {
79d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr);
80d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
81d132466eSBarry Smith         ierr = ViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
829b94acceSBarry Smith       }
839b94acceSBarry Smith     }
840f5bd95cSBarry Smith   } else if (isstring) {
85454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
86454a90a3SBarry Smith     ierr = ViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
8719bcc07fSBarry Smith   }
8877ed5343SBarry Smith   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
890ef38995SBarry Smith   ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr);
909b94acceSBarry Smith   ierr = SLESView(sles,viewer);CHKERRQ(ierr);
910ef38995SBarry Smith   ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr);
923a40ed3dSBarry Smith   PetscFunctionReturn(0);
939b94acceSBarry Smith }
949b94acceSBarry Smith 
9515091d37SBarry Smith #undef __FUNC__
96b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFromOptions"
979b94acceSBarry Smith /*@
98682d7d0cSBarry Smith    SNESSetFromOptions - Sets various SNES and SLES parameters from user options.
999b94acceSBarry Smith 
100c7afd0dbSLois Curfman McInnes    Collective on SNES
101c7afd0dbSLois Curfman McInnes 
1029b94acceSBarry Smith    Input Parameter:
1039b94acceSBarry Smith .  snes - the SNES context
1049b94acceSBarry Smith 
10536851e7fSLois Curfman McInnes    Options Database Keys:
1066831982aSBarry Smith +  -snes_type <type> - ls, tr, umls, umtr, test
10782738288SBarry Smith .  -snes_stol - convergence tolerance in terms of the norm
10882738288SBarry Smith                 of the change in the solution between steps
109b39c3a46SLois Curfman McInnes .  -snes_atol <atol> - absolute tolerance of residual norm
110b39c3a46SLois Curfman McInnes .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
111b39c3a46SLois Curfman McInnes .  -snes_max_it <max_it> - maximum number of iterations
112b39c3a46SLois Curfman McInnes .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
113b39c3a46SLois Curfman McInnes .  -snes_trtol <trtol> - trust region tolerance
11482738288SBarry Smith .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
11582738288SBarry Smith                                solver; hence iterations will continue until max_it
1161fbbfb26SLois Curfman McInnes                                or some other criterion is reached. Saves expense
11782738288SBarry Smith                                of convergence test
11882738288SBarry Smith .  -snes_monitor - prints residual norm at each iteration
119d132466eSBarry Smith .  -snes_vecmonitor - plots solution at each iteration
120d132466eSBarry Smith .  -snes_vecmonitor_update - plots update to solution at each iteration
12182738288SBarry Smith .  -snes_xmonitor - plots residual norm at each iteration
122e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
12336851e7fSLois Curfman McInnes -  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
12482738288SBarry Smith 
12582738288SBarry Smith     Options Database for Eisenstat-Walker method:
12682738288SBarry Smith +  -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence
12782738288SBarry Smith .  -snes_ksp_eq_version ver - version of  Eisenstat-Walker method
12836851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
12936851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
13036851e7fSLois Curfman McInnes .  -snes_ksp_ew_gamma <gamma> - Sets gamma
13136851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha <alpha> - Sets alpha
13236851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
13336851e7fSLois Curfman McInnes -  -snes_ksp_ew_threshold <threshold> - Sets threshold
13482738288SBarry Smith 
13511ca99fdSLois Curfman McInnes    Notes:
13611ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
13711ca99fdSLois Curfman McInnes    the users manual.
13883e2fdc7SBarry Smith 
13936851e7fSLois Curfman McInnes    Level: beginner
14036851e7fSLois Curfman McInnes 
1419b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1429b94acceSBarry Smith 
143186905e3SBarry Smith .seealso: SNESSetOptionsPrefix(), SNESSetTypeFromOptions()
1449b94acceSBarry Smith @*/
1459b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1469b94acceSBarry Smith {
1479b94acceSBarry Smith   SLES                sles;
148186905e3SBarry Smith   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
149f1af5d2fSBarry Smith   PetscTruth          flg;
150186905e3SBarry Smith   int                 ierr;
151186905e3SBarry Smith   char                *deft,type[256];
1529b94acceSBarry Smith 
1533a40ed3dSBarry Smith   PetscFunctionBegin;
15477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
155ca161407SBarry Smith 
156cd5578b5SBarry Smith   ierr = OptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
157186905e3SBarry Smith     if (snes->type_name) {
158186905e3SBarry Smith       deft = snes->type_name;
159186905e3SBarry Smith     } else {
160186905e3SBarry Smith       if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
161186905e3SBarry Smith         deft = SNESEQLS;
162186905e3SBarry Smith       } else {
163186905e3SBarry Smith         deft = SNESUMTR;
164d64ed03dSBarry Smith       }
165d64ed03dSBarry Smith     }
1664bbc92c1SBarry Smith 
167186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
168186905e3SBarry Smith     ierr = OptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
169d64ed03dSBarry Smith     if (flg) {
170186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
171186905e3SBarry Smith     } else if (!snes->type_name) {
172186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
173d64ed03dSBarry Smith     }
17493c39befSBarry Smith 
175186905e3SBarry Smith     ierr = OptionsDouble("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
176186905e3SBarry Smith     ierr = OptionsDouble("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr);
177186905e3SBarry Smith 
178186905e3SBarry Smith     ierr = OptionsDouble("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
179186905e3SBarry Smith     ierr = OptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
180186905e3SBarry Smith     ierr = OptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
181186905e3SBarry Smith     ierr = OptionsDouble("-snes_fmin","Minimization function tolerance","SNESSetMinimizationFunctionTolerance",snes->fmin,&snes->fmin,0);CHKERRQ(ierr);
182186905e3SBarry Smith 
183186905e3SBarry Smith     ierr = OptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr);
184186905e3SBarry Smith 
185186905e3SBarry Smith     ierr = OptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
186186905e3SBarry Smith     ierr = OptionsDouble("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
187186905e3SBarry Smith     ierr = OptionsDouble("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
188186905e3SBarry Smith     ierr = OptionsDouble("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
189186905e3SBarry Smith     ierr = OptionsDouble("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
190186905e3SBarry Smith     ierr = OptionsDouble("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
191186905e3SBarry Smith     ierr = OptionsDouble("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
192186905e3SBarry Smith 
193186905e3SBarry Smith     ierr = OptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
19493c39befSBarry Smith     if (flg) {snes->converged = 0;}
195186905e3SBarry Smith     ierr = OptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
1965cd90555SBarry Smith     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
197186905e3SBarry Smith     ierr = OptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
198b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
199186905e3SBarry Smith     ierr = OptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
200b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
201186905e3SBarry Smith     ierr = OptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
202b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
203186905e3SBarry Smith     ierr = OptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
2047c922b88SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
205186905e3SBarry Smith     ierr = OptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
206186905e3SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
207e24b481bSBarry Smith 
208186905e3SBarry Smith     ierr = OptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
2093c7409f5SSatish Balay     if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) {
210186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
211981c4779SBarry Smith       PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
212981c4779SBarry Smith     } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
213186905e3SBarry Smith       ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeHessian,snes->funP);CHKERRQ(ierr);
214d64ed03dSBarry Smith       PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n");
2159b94acceSBarry Smith     }
216639f9d9dSBarry Smith 
217186905e3SBarry Smith     if (snes->setfromoptions) {
218186905e3SBarry Smith       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
219639f9d9dSBarry Smith     }
220639f9d9dSBarry Smith 
2214bbc92c1SBarry Smith   ierr = OptionsEnd();CHKERRQ(ierr);
2224bbc92c1SBarry Smith 
2239b94acceSBarry Smith   ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
2249b94acceSBarry Smith   ierr = SLESSetFromOptions(sles);CHKERRQ(ierr);
22593993e2dSLois Curfman McInnes 
2263a40ed3dSBarry Smith   PetscFunctionReturn(0);
2279b94acceSBarry Smith }
2289b94acceSBarry Smith 
229a847f771SSatish Balay 
2305615d1e5SSatish Balay #undef __FUNC__
231b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetApplicationContext"
2329b94acceSBarry Smith /*@
2339b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2349b94acceSBarry Smith    the nonlinear solvers.
2359b94acceSBarry Smith 
236fee21e36SBarry Smith    Collective on SNES
237fee21e36SBarry Smith 
238c7afd0dbSLois Curfman McInnes    Input Parameters:
239c7afd0dbSLois Curfman McInnes +  snes - the SNES context
240c7afd0dbSLois Curfman McInnes -  usrP - optional user context
241c7afd0dbSLois Curfman McInnes 
24236851e7fSLois Curfman McInnes    Level: intermediate
24336851e7fSLois Curfman McInnes 
2449b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2459b94acceSBarry Smith 
2469b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2479b94acceSBarry Smith @*/
2489b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
2499b94acceSBarry Smith {
2503a40ed3dSBarry Smith   PetscFunctionBegin;
25177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2529b94acceSBarry Smith   snes->user		= usrP;
2533a40ed3dSBarry Smith   PetscFunctionReturn(0);
2549b94acceSBarry Smith }
25574679c65SBarry Smith 
2565615d1e5SSatish Balay #undef __FUNC__
257b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetApplicationContext"
2589b94acceSBarry Smith /*@C
2599b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
2609b94acceSBarry Smith    nonlinear solvers.
2619b94acceSBarry Smith 
262c7afd0dbSLois Curfman McInnes    Not Collective
263c7afd0dbSLois Curfman McInnes 
2649b94acceSBarry Smith    Input Parameter:
2659b94acceSBarry Smith .  snes - SNES context
2669b94acceSBarry Smith 
2679b94acceSBarry Smith    Output Parameter:
2689b94acceSBarry Smith .  usrP - user context
2699b94acceSBarry Smith 
27036851e7fSLois Curfman McInnes    Level: intermediate
27136851e7fSLois Curfman McInnes 
2729b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
2739b94acceSBarry Smith 
2749b94acceSBarry Smith .seealso: SNESSetApplicationContext()
2759b94acceSBarry Smith @*/
2769b94acceSBarry Smith int SNESGetApplicationContext(SNES snes,void **usrP)
2779b94acceSBarry Smith {
2783a40ed3dSBarry Smith   PetscFunctionBegin;
27977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2809b94acceSBarry Smith   *usrP = snes->user;
2813a40ed3dSBarry Smith   PetscFunctionReturn(0);
2829b94acceSBarry Smith }
28374679c65SBarry Smith 
2845615d1e5SSatish Balay #undef __FUNC__
285b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetIterationNumber"
2869b94acceSBarry Smith /*@
287c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
288c8228a4eSBarry Smith    at this time.
2899b94acceSBarry Smith 
290c7afd0dbSLois Curfman McInnes    Not Collective
291c7afd0dbSLois Curfman McInnes 
2929b94acceSBarry Smith    Input Parameter:
2939b94acceSBarry Smith .  snes - SNES context
2949b94acceSBarry Smith 
2959b94acceSBarry Smith    Output Parameter:
2969b94acceSBarry Smith .  iter - iteration number
2979b94acceSBarry Smith 
298c8228a4eSBarry Smith    Notes:
299c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
300c8228a4eSBarry Smith 
301c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
30208405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
30308405cd6SLois Curfman McInnes .vb
30408405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
30508405cd6SLois Curfman McInnes       if (!(it % 2)) {
30608405cd6SLois Curfman McInnes         [compute Jacobian here]
30708405cd6SLois Curfman McInnes       }
30808405cd6SLois Curfman McInnes .ve
309c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
31008405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
311c8228a4eSBarry Smith 
31236851e7fSLois Curfman McInnes    Level: intermediate
31336851e7fSLois Curfman McInnes 
3149b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3159b94acceSBarry Smith @*/
3169b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3179b94acceSBarry Smith {
3183a40ed3dSBarry Smith   PetscFunctionBegin;
31977c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
32074679c65SBarry Smith   PetscValidIntPointer(iter);
3219b94acceSBarry Smith   *iter = snes->iter;
3223a40ed3dSBarry Smith   PetscFunctionReturn(0);
3239b94acceSBarry Smith }
32474679c65SBarry Smith 
3255615d1e5SSatish Balay #undef __FUNC__
326b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetFunctionNorm"
3279b94acceSBarry Smith /*@
3289b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3299b94acceSBarry Smith    with SNESSSetFunction().
3309b94acceSBarry Smith 
331c7afd0dbSLois Curfman McInnes    Collective on SNES
332c7afd0dbSLois Curfman McInnes 
3339b94acceSBarry Smith    Input Parameter:
3349b94acceSBarry Smith .  snes - SNES context
3359b94acceSBarry Smith 
3369b94acceSBarry Smith    Output Parameter:
3379b94acceSBarry Smith .  fnorm - 2-norm of function
3389b94acceSBarry Smith 
3399b94acceSBarry Smith    Note:
3409b94acceSBarry Smith    SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only.
341a86d99e1SLois Curfman McInnes    A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is
342a86d99e1SLois Curfman McInnes    SNESGetGradientNorm().
3439b94acceSBarry Smith 
34436851e7fSLois Curfman McInnes    Level: intermediate
34536851e7fSLois Curfman McInnes 
3469b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
347a86d99e1SLois Curfman McInnes 
34808405cd6SLois Curfman McInnes .seealso: SNESGetFunction()
3499b94acceSBarry Smith @*/
3509b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm)
3519b94acceSBarry Smith {
3523a40ed3dSBarry Smith   PetscFunctionBegin;
35377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
35474679c65SBarry Smith   PetscValidScalarPointer(fnorm);
35574679c65SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
356*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_NONLINEAR_EQUATIONS only");
35774679c65SBarry Smith   }
3589b94acceSBarry Smith   *fnorm = snes->norm;
3593a40ed3dSBarry Smith   PetscFunctionReturn(0);
3609b94acceSBarry Smith }
36174679c65SBarry Smith 
3625615d1e5SSatish Balay #undef __FUNC__
363b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetGradientNorm"
3649b94acceSBarry Smith /*@
3659b94acceSBarry Smith    SNESGetGradientNorm - Gets the norm of the current gradient that was set
3669b94acceSBarry Smith    with SNESSSetGradient().
3679b94acceSBarry Smith 
368c7afd0dbSLois Curfman McInnes    Collective on SNES
369c7afd0dbSLois Curfman McInnes 
3709b94acceSBarry Smith    Input Parameter:
3719b94acceSBarry Smith .  snes - SNES context
3729b94acceSBarry Smith 
3739b94acceSBarry Smith    Output Parameter:
3749b94acceSBarry Smith .  fnorm - 2-norm of gradient
3759b94acceSBarry Smith 
3769b94acceSBarry Smith    Note:
3779b94acceSBarry Smith    SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION
378a86d99e1SLois Curfman McInnes    methods only.  A related routine for SNES_NONLINEAR_EQUATIONS methods
379a86d99e1SLois Curfman McInnes    is SNESGetFunctionNorm().
3809b94acceSBarry Smith 
38136851e7fSLois Curfman McInnes    Level: intermediate
38236851e7fSLois Curfman McInnes 
3839b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm
384a86d99e1SLois Curfman McInnes 
385a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient()
3869b94acceSBarry Smith @*/
3879b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm)
3889b94acceSBarry Smith {
3893a40ed3dSBarry Smith   PetscFunctionBegin;
39077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
39174679c65SBarry Smith   PetscValidScalarPointer(gnorm);
39274679c65SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
393*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"For SNES_UNCONSTRAINED_MINIMIZATION only");
39474679c65SBarry Smith   }
3959b94acceSBarry Smith   *gnorm = snes->norm;
3963a40ed3dSBarry Smith   PetscFunctionReturn(0);
3979b94acceSBarry Smith }
39874679c65SBarry Smith 
3995615d1e5SSatish Balay #undef __FUNC__
400b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetNumberUnsuccessfulSteps"
4019b94acceSBarry Smith /*@
4029b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4039b94acceSBarry Smith    attempted by the nonlinear solver.
4049b94acceSBarry Smith 
405c7afd0dbSLois Curfman McInnes    Not Collective
406c7afd0dbSLois Curfman McInnes 
4079b94acceSBarry Smith    Input Parameter:
4089b94acceSBarry Smith .  snes - SNES context
4099b94acceSBarry Smith 
4109b94acceSBarry Smith    Output Parameter:
4119b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4129b94acceSBarry Smith 
413c96a6f78SLois Curfman McInnes    Notes:
414c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
415c96a6f78SLois Curfman McInnes 
41636851e7fSLois Curfman McInnes    Level: intermediate
41736851e7fSLois Curfman McInnes 
4189b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4199b94acceSBarry Smith @*/
4209b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4219b94acceSBarry Smith {
4223a40ed3dSBarry Smith   PetscFunctionBegin;
42377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
42474679c65SBarry Smith   PetscValidIntPointer(nfails);
4259b94acceSBarry Smith   *nfails = snes->nfailures;
4263a40ed3dSBarry Smith   PetscFunctionReturn(0);
4279b94acceSBarry Smith }
428a847f771SSatish Balay 
4295615d1e5SSatish Balay #undef __FUNC__
430b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetNumberLinearIterations"
431c96a6f78SLois Curfman McInnes /*@
432c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
433c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
434c96a6f78SLois Curfman McInnes 
435c7afd0dbSLois Curfman McInnes    Not Collective
436c7afd0dbSLois Curfman McInnes 
437c96a6f78SLois Curfman McInnes    Input Parameter:
438c96a6f78SLois Curfman McInnes .  snes - SNES context
439c96a6f78SLois Curfman McInnes 
440c96a6f78SLois Curfman McInnes    Output Parameter:
441c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
442c96a6f78SLois Curfman McInnes 
443c96a6f78SLois Curfman McInnes    Notes:
444c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
445c96a6f78SLois Curfman McInnes 
44636851e7fSLois Curfman McInnes    Level: intermediate
44736851e7fSLois Curfman McInnes 
448c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
449c96a6f78SLois Curfman McInnes @*/
45086bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
451c96a6f78SLois Curfman McInnes {
4523a40ed3dSBarry Smith   PetscFunctionBegin;
453c96a6f78SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
454c96a6f78SLois Curfman McInnes   PetscValidIntPointer(lits);
455c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
4563a40ed3dSBarry Smith   PetscFunctionReturn(0);
457c96a6f78SLois Curfman McInnes }
458c96a6f78SLois Curfman McInnes 
459c96a6f78SLois Curfman McInnes #undef __FUNC__
460b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetSLES"
4619b94acceSBarry Smith /*@C
4629b94acceSBarry Smith    SNESGetSLES - Returns the SLES context for a SNES solver.
4639b94acceSBarry Smith 
464c7afd0dbSLois Curfman McInnes    Not Collective, but if SNES object is parallel, then SLES object is parallel
465c7afd0dbSLois Curfman McInnes 
4669b94acceSBarry Smith    Input Parameter:
4679b94acceSBarry Smith .  snes - the SNES context
4689b94acceSBarry Smith 
4699b94acceSBarry Smith    Output Parameter:
4709b94acceSBarry Smith .  sles - the SLES context
4719b94acceSBarry Smith 
4729b94acceSBarry Smith    Notes:
4739b94acceSBarry Smith    The user can then directly manipulate the SLES context to set various
4749b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
4759b94acceSBarry Smith    KSP and PC contexts as well.
4769b94acceSBarry Smith 
47736851e7fSLois Curfman McInnes    Level: beginner
47836851e7fSLois Curfman McInnes 
4799b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context
4809b94acceSBarry Smith 
4819b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP()
4829b94acceSBarry Smith @*/
4839b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles)
4849b94acceSBarry Smith {
4853a40ed3dSBarry Smith   PetscFunctionBegin;
48677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
4879b94acceSBarry Smith   *sles = snes->sles;
4883a40ed3dSBarry Smith   PetscFunctionReturn(0);
4899b94acceSBarry Smith }
49082bf6240SBarry Smith 
491e24b481bSBarry Smith #undef __FUNC__
492b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESPublish_Petsc"
493454a90a3SBarry Smith static int SNESPublish_Petsc(PetscObject obj)
494e24b481bSBarry Smith {
495aa482453SBarry Smith #if defined(PETSC_HAVE_AMS)
496454a90a3SBarry Smith   SNES          v = (SNES) obj;
497e24b481bSBarry Smith   int          ierr;
49843d6d2cbSBarry Smith #endif
499e24b481bSBarry Smith 
500e24b481bSBarry Smith   PetscFunctionBegin;
501e24b481bSBarry Smith 
50243d6d2cbSBarry Smith #if defined(PETSC_HAVE_AMS)
503e24b481bSBarry Smith   /* if it is already published then return */
504e24b481bSBarry Smith   if (v->amem >=0) PetscFunctionReturn(0);
505e24b481bSBarry Smith 
506454a90a3SBarry Smith   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
507e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
508e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
509e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
510e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
511454a90a3SBarry Smith   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
512433994e6SBarry Smith #endif
513e24b481bSBarry Smith   PetscFunctionReturn(0);
514e24b481bSBarry Smith }
515e24b481bSBarry Smith 
5169b94acceSBarry Smith /* -----------------------------------------------------------*/
5175615d1e5SSatish Balay #undef __FUNC__
518b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESCreate"
5199b94acceSBarry Smith /*@C
5209b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5219b94acceSBarry Smith 
522c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
523c7afd0dbSLois Curfman McInnes 
524c7afd0dbSLois Curfman McInnes    Input Parameters:
525c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
526c7afd0dbSLois Curfman McInnes -  type - type of method, either
527c7afd0dbSLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations)
528c7afd0dbSLois Curfman McInnes    or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization)
5299b94acceSBarry Smith 
5309b94acceSBarry Smith    Output Parameter:
5319b94acceSBarry Smith .  outsnes - the new SNES context
5329b94acceSBarry Smith 
533c7afd0dbSLois Curfman McInnes    Options Database Keys:
534c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
535c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
536c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
537c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
538c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
539c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
540c1f60f51SBarry Smith 
54136851e7fSLois Curfman McInnes    Level: beginner
54236851e7fSLois Curfman McInnes 
5439b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5449b94acceSBarry Smith 
54563a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy()
5469b94acceSBarry Smith @*/
5474b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes)
5489b94acceSBarry Smith {
5499b94acceSBarry Smith   int                 ierr;
5509b94acceSBarry Smith   SNES                snes;
5519b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
55237fcc0dbSBarry Smith 
5533a40ed3dSBarry Smith   PetscFunctionBegin;
5549b94acceSBarry Smith   *outsnes = 0;
555d64ed03dSBarry Smith   if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){
556*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"incorrect method type");
557d64ed03dSBarry Smith   }
5583f1db9ecSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);
5599b94acceSBarry Smith   PLogObjectCreate(snes);
560e24b481bSBarry Smith   snes->bops->publish     = SNESPublish_Petsc;
5619b94acceSBarry Smith   snes->max_its           = 50;
5629750a799SBarry Smith   snes->max_funcs	  = 10000;
5639b94acceSBarry Smith   snes->norm		  = 0.0;
5645a2d0531SBarry Smith   if (type == SNES_UNCONSTRAINED_MINIMIZATION) {
565b18e04deSLois Curfman McInnes     snes->rtol		  = 1.e-8;
566b18e04deSLois Curfman McInnes     snes->ttol            = 0.0;
5679b94acceSBarry Smith     snes->atol		  = 1.e-10;
568d64ed03dSBarry Smith   } else {
569b4874afaSBarry Smith     snes->rtol		  = 1.e-8;
570b4874afaSBarry Smith     snes->ttol            = 0.0;
571b4874afaSBarry Smith     snes->atol		  = 1.e-50;
572b4874afaSBarry Smith   }
5739b94acceSBarry Smith   snes->xtol		  = 1.e-8;
574b18e04deSLois Curfman McInnes   snes->trunctol	  = 1.e-12; /* no longer used */
5759b94acceSBarry Smith   snes->nfuncs            = 0;
5769b94acceSBarry Smith   snes->nfailures         = 0;
5777a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
578639f9d9dSBarry Smith   snes->numbermonitors    = 0;
5799b94acceSBarry Smith   snes->data              = 0;
5809b94acceSBarry Smith   snes->view              = 0;
5819b94acceSBarry Smith   snes->computeumfunction = 0;
5829b94acceSBarry Smith   snes->umfunP            = 0;
5839b94acceSBarry Smith   snes->fc                = 0;
5849b94acceSBarry Smith   snes->deltatol          = 1.e-12;
5859b94acceSBarry Smith   snes->fmin              = -1.e30;
5869b94acceSBarry Smith   snes->method_class      = type;
5879b94acceSBarry Smith   snes->set_method_called = 0;
58882bf6240SBarry Smith   snes->setupcalled       = 0;
589186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
5906f24a144SLois Curfman McInnes   snes->vwork             = 0;
5916f24a144SLois Curfman McInnes   snes->nwork             = 0;
592758f92a0SBarry Smith   snes->conv_hist_len     = 0;
593758f92a0SBarry Smith   snes->conv_hist_max     = 0;
594758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
595758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
596758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
597184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
5989b94acceSBarry Smith 
5999b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
6000452661fSBarry Smith   kctx = PetscNew(SNES_KSP_EW_ConvCtx);CHKPTRQ(kctx);
601eed86810SBarry Smith   PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
6029b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6039b94acceSBarry Smith   kctx->version     = 2;
6049b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6059b94acceSBarry Smith                              this was too large for some test cases */
6069b94acceSBarry Smith   kctx->rtol_last   = 0;
6079b94acceSBarry Smith   kctx->rtol_max    = .9;
6089b94acceSBarry Smith   kctx->gamma       = 1.0;
6099b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6109b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6119b94acceSBarry Smith   kctx->threshold   = .1;
6129b94acceSBarry Smith   kctx->lresid_last = 0;
6139b94acceSBarry Smith   kctx->norm_last   = 0;
6149b94acceSBarry Smith 
6159b94acceSBarry Smith   ierr = SLESCreate(comm,&snes->sles);CHKERRQ(ierr);
6169b94acceSBarry Smith   PLogObjectParent(snes,snes->sles)
6175334005bSBarry Smith 
6189b94acceSBarry Smith   *outsnes = snes;
61900036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
6203a40ed3dSBarry Smith   PetscFunctionReturn(0);
6219b94acceSBarry Smith }
6229b94acceSBarry Smith 
6239b94acceSBarry Smith /* --------------------------------------------------------------- */
6245615d1e5SSatish Balay #undef __FUNC__
625b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetFunction"
6269b94acceSBarry Smith /*@C
6279b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6289b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6299b94acceSBarry Smith    equations.
6309b94acceSBarry Smith 
631fee21e36SBarry Smith    Collective on SNES
632fee21e36SBarry Smith 
633c7afd0dbSLois Curfman McInnes    Input Parameters:
634c7afd0dbSLois Curfman McInnes +  snes - the SNES context
635c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
636c7afd0dbSLois Curfman McInnes .  r - vector to store function value
637c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
638c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6399b94acceSBarry Smith 
640c7afd0dbSLois Curfman McInnes    Calling sequence of func:
6418d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
642c7afd0dbSLois Curfman McInnes 
643313e4042SLois Curfman McInnes .  f - function vector
644c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
6459b94acceSBarry Smith 
6469b94acceSBarry Smith    Notes:
6479b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6489b94acceSBarry Smith $      f'(x) x = -f(x),
649c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
6509b94acceSBarry Smith 
6519b94acceSBarry Smith    SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6529b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6539b94acceSBarry Smith    SNESSetMinimizationFunction() and SNESSetGradient();
6549b94acceSBarry Smith 
65536851e7fSLois Curfman McInnes    Level: beginner
65636851e7fSLois Curfman McInnes 
6579b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6589b94acceSBarry Smith 
659a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6609b94acceSBarry Smith @*/
6615334005bSBarry Smith int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
6629b94acceSBarry Smith {
6633a40ed3dSBarry Smith   PetscFunctionBegin;
66477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
665184914b5SBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE);
666184914b5SBarry Smith   PetscCheckSameComm(snes,r);
667a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
668*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
669a8c6a408SBarry Smith   }
670184914b5SBarry Smith 
6719b94acceSBarry Smith   snes->computefunction     = func;
6729b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
6739b94acceSBarry Smith   snes->funP                = ctx;
6743a40ed3dSBarry Smith   PetscFunctionReturn(0);
6759b94acceSBarry Smith }
6769b94acceSBarry Smith 
6775615d1e5SSatish Balay #undef __FUNC__
678b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESComputeFunction"
6799b94acceSBarry Smith /*@
68036851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
6819b94acceSBarry Smith                          SNESSetFunction().
6829b94acceSBarry Smith 
683c7afd0dbSLois Curfman McInnes    Collective on SNES
684c7afd0dbSLois Curfman McInnes 
6859b94acceSBarry Smith    Input Parameters:
686c7afd0dbSLois Curfman McInnes +  snes - the SNES context
687c7afd0dbSLois Curfman McInnes -  x - input vector
6889b94acceSBarry Smith 
6899b94acceSBarry Smith    Output Parameter:
6903638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
6919b94acceSBarry Smith 
6921bffabb2SLois Curfman McInnes    Notes:
6939b94acceSBarry Smith    SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only.
6949b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
6959b94acceSBarry Smith    SNESComputeMinimizationFunction() and SNESComputeGradient();
6969b94acceSBarry Smith 
69736851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
69836851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
69936851e7fSLois Curfman McInnes    themselves.
70036851e7fSLois Curfman McInnes 
70136851e7fSLois Curfman McInnes    Level: developer
70236851e7fSLois Curfman McInnes 
7039b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7049b94acceSBarry Smith 
705a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7069b94acceSBarry Smith @*/
7079b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x,Vec y)
7089b94acceSBarry Smith {
7099b94acceSBarry Smith   int    ierr;
7109b94acceSBarry Smith 
7113a40ed3dSBarry Smith   PetscFunctionBegin;
712184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
713184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
714184914b5SBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE);
715184914b5SBarry Smith   PetscCheckSameComm(snes,x);
716184914b5SBarry Smith   PetscCheckSameComm(snes,y);
717d4bb536fSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
718*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
719d4bb536fSBarry Smith   }
720184914b5SBarry Smith 
72181bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
722d64ed03dSBarry Smith   PetscStackPush("SNES user function");
7239b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
724d64ed03dSBarry Smith   PetscStackPop;
725ae3c334cSLois Curfman McInnes   snes->nfuncs++;
72681bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
7273a40ed3dSBarry Smith   PetscFunctionReturn(0);
7289b94acceSBarry Smith }
7299b94acceSBarry Smith 
7305615d1e5SSatish Balay #undef __FUNC__
731b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetMinimizationFunction"
7329b94acceSBarry Smith /*@C
7339b94acceSBarry Smith    SNESSetMinimizationFunction - Sets the function evaluation routine for
7349b94acceSBarry Smith    unconstrained minimization.
7359b94acceSBarry Smith 
736fee21e36SBarry Smith    Collective on SNES
737fee21e36SBarry Smith 
738c7afd0dbSLois Curfman McInnes    Input Parameters:
739c7afd0dbSLois Curfman McInnes +  snes - the SNES context
740c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
741c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
742c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
7439b94acceSBarry Smith 
744c7afd0dbSLois Curfman McInnes    Calling sequence of func:
745329f5518SBarry Smith $     func (SNES snes,Vec x,PetscReal *f,void *ctx);
746c7afd0dbSLois Curfman McInnes 
747c7afd0dbSLois Curfman McInnes +  x - input vector
7489b94acceSBarry Smith .  f - function
749c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined function context
7509b94acceSBarry Smith 
75136851e7fSLois Curfman McInnes    Level: beginner
75236851e7fSLois Curfman McInnes 
7539b94acceSBarry Smith    Notes:
7549b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
7559b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
7569b94acceSBarry Smith    SNESSetFunction().
7579b94acceSBarry Smith 
7589b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function
7599b94acceSBarry Smith 
760a86d99e1SLois Curfman McInnes .seealso:  SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(),
761a86d99e1SLois Curfman McInnes            SNESSetHessian(), SNESSetGradient()
7629b94acceSBarry Smith @*/
763329f5518SBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,PetscReal*,void*),void *ctx)
7649b94acceSBarry Smith {
7653a40ed3dSBarry Smith   PetscFunctionBegin;
76677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
767a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
768*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
769a8c6a408SBarry Smith   }
7709b94acceSBarry Smith   snes->computeumfunction   = func;
7719b94acceSBarry Smith   snes->umfunP              = ctx;
7723a40ed3dSBarry Smith   PetscFunctionReturn(0);
7739b94acceSBarry Smith }
7749b94acceSBarry Smith 
7755615d1e5SSatish Balay #undef __FUNC__
776b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESComputeMinimizationFunction"
7779b94acceSBarry Smith /*@
7789b94acceSBarry Smith    SNESComputeMinimizationFunction - Computes the function that has been
7799b94acceSBarry Smith    set with SNESSetMinimizationFunction().
7809b94acceSBarry Smith 
781c7afd0dbSLois Curfman McInnes    Collective on SNES
782c7afd0dbSLois Curfman McInnes 
7839b94acceSBarry Smith    Input Parameters:
784c7afd0dbSLois Curfman McInnes +  snes - the SNES context
785c7afd0dbSLois Curfman McInnes -  x - input vector
7869b94acceSBarry Smith 
7879b94acceSBarry Smith    Output Parameter:
7889b94acceSBarry Smith .  y - function value
7899b94acceSBarry Smith 
7909b94acceSBarry Smith    Notes:
7919b94acceSBarry Smith    SNESComputeMinimizationFunction() is valid only for
7929b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
7939b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
794a86d99e1SLois Curfman McInnes 
79536851e7fSLois Curfman McInnes    SNESComputeMinimizationFunction() is typically used within minimization
79636851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
79736851e7fSLois Curfman McInnes    themselves.
79836851e7fSLois Curfman McInnes 
79936851e7fSLois Curfman McInnes    Level: developer
80036851e7fSLois Curfman McInnes 
801a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function
802a86d99e1SLois Curfman McInnes 
803a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(),
804a86d99e1SLois Curfman McInnes           SNESComputeGradient(), SNESComputeHessian()
8059b94acceSBarry Smith @*/
806329f5518SBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,PetscReal *y)
8079b94acceSBarry Smith {
8089b94acceSBarry Smith   int    ierr;
8093a40ed3dSBarry Smith 
8103a40ed3dSBarry Smith   PetscFunctionBegin;
811184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
812184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
813184914b5SBarry Smith   PetscCheckSameComm(snes,x);
814a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
815*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"Only for SNES_UNCONSTRAINED_MINIMIZATION");
816a8c6a408SBarry Smith   }
817184914b5SBarry Smith 
81881bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr);
819d64ed03dSBarry Smith   PetscStackPush("SNES user minimzation function");
8209b94acceSBarry Smith   ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP);CHKERRQ(ierr);
821d64ed03dSBarry Smith   PetscStackPop;
822ae3c334cSLois Curfman McInnes   snes->nfuncs++;
82381bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0);CHKERRQ(ierr);
8243a40ed3dSBarry Smith   PetscFunctionReturn(0);
8259b94acceSBarry Smith }
8269b94acceSBarry Smith 
8275615d1e5SSatish Balay #undef __FUNC__
828b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetGradient"
8299b94acceSBarry Smith /*@C
8309b94acceSBarry Smith    SNESSetGradient - Sets the gradient evaluation routine and gradient
8319b94acceSBarry Smith    vector for use by the SNES routines.
8329b94acceSBarry Smith 
833c7afd0dbSLois Curfman McInnes    Collective on SNES
834c7afd0dbSLois Curfman McInnes 
8359b94acceSBarry Smith    Input Parameters:
836c7afd0dbSLois Curfman McInnes +  snes - the SNES context
8379b94acceSBarry Smith .  func - function evaluation routine
838044dda88SLois Curfman McInnes .  ctx - optional user-defined context for private data for the
839044dda88SLois Curfman McInnes          gradient evaluation routine (may be PETSC_NULL)
840c7afd0dbSLois Curfman McInnes -  r - vector to store gradient value
841fee21e36SBarry Smith 
8429b94acceSBarry Smith    Calling sequence of func:
8438d76a1e5SLois Curfman McInnes $     func (SNES, Vec x, Vec g, void *ctx);
8449b94acceSBarry Smith 
845c7afd0dbSLois Curfman McInnes +  x - input vector
8469b94acceSBarry Smith .  g - gradient vector
847c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined gradient context
8489b94acceSBarry Smith 
8499b94acceSBarry Smith    Notes:
8509b94acceSBarry Smith    SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
8519b94acceSBarry Smith    methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
8529b94acceSBarry Smith    SNESSetFunction().
8539b94acceSBarry Smith 
85436851e7fSLois Curfman McInnes    Level: beginner
85536851e7fSLois Curfman McInnes 
8569b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
8579b94acceSBarry Smith 
858a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(),
859a86d99e1SLois Curfman McInnes           SNESSetMinimizationFunction(),
8609b94acceSBarry Smith @*/
86174679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
8629b94acceSBarry Smith {
8633a40ed3dSBarry Smith   PetscFunctionBegin;
86477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
865184914b5SBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE);
866184914b5SBarry Smith   PetscCheckSameComm(snes,r);
867a8c6a408SBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
868*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
869a8c6a408SBarry Smith   }
8709b94acceSBarry Smith   snes->computefunction     = func;
8719b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
8729b94acceSBarry Smith   snes->funP                = ctx;
8733a40ed3dSBarry Smith   PetscFunctionReturn(0);
8749b94acceSBarry Smith }
8759b94acceSBarry Smith 
8765615d1e5SSatish Balay #undef __FUNC__
877b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESComputeGradient"
8789b94acceSBarry Smith /*@
879a86d99e1SLois Curfman McInnes    SNESComputeGradient - Computes the gradient that has been set with
880a86d99e1SLois Curfman McInnes    SNESSetGradient().
8819b94acceSBarry Smith 
882c7afd0dbSLois Curfman McInnes    Collective on SNES
883c7afd0dbSLois Curfman McInnes 
8849b94acceSBarry Smith    Input Parameters:
885c7afd0dbSLois Curfman McInnes +  snes - the SNES context
886c7afd0dbSLois Curfman McInnes -  x - input vector
8879b94acceSBarry Smith 
8889b94acceSBarry Smith    Output Parameter:
8899b94acceSBarry Smith .  y - gradient vector
8909b94acceSBarry Smith 
8919b94acceSBarry Smith    Notes:
8929b94acceSBarry Smith    SNESComputeGradient() is valid only for
8939b94acceSBarry Smith    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
8949b94acceSBarry Smith    SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction().
895a86d99e1SLois Curfman McInnes 
89636851e7fSLois Curfman McInnes    SNESComputeGradient() is typically used within minimization
89736851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
89836851e7fSLois Curfman McInnes    themselves.
89936851e7fSLois Curfman McInnes 
90036851e7fSLois Curfman McInnes    Level: developer
90136851e7fSLois Curfman McInnes 
902a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient
903a86d99e1SLois Curfman McInnes 
904a86d99e1SLois Curfman McInnes .seealso:  SNESSetGradient(), SNESGetGradient(),
905a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction(), SNESComputeHessian()
9069b94acceSBarry Smith @*/
9079b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x,Vec y)
9089b94acceSBarry Smith {
9099b94acceSBarry Smith   int    ierr;
9103a40ed3dSBarry Smith 
9113a40ed3dSBarry Smith   PetscFunctionBegin;
912184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
913184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
914184914b5SBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE);
915184914b5SBarry Smith   PetscCheckSameComm(snes,x);
916184914b5SBarry Smith   PetscCheckSameComm(snes,y);
9173a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
918*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
9193a40ed3dSBarry Smith   }
92081bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr);
921d64ed03dSBarry Smith   PetscStackPush("SNES user gradient function");
9229b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
923d64ed03dSBarry Smith   PetscStackPop;
92481bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_GradientEval,snes,x,y,0);CHKERRQ(ierr);
9253a40ed3dSBarry Smith   PetscFunctionReturn(0);
9269b94acceSBarry Smith }
9279b94acceSBarry Smith 
9285615d1e5SSatish Balay #undef __FUNC__
929b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESComputeJacobian"
93062fef451SLois Curfman McInnes /*@
93162fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
93262fef451SLois Curfman McInnes    set with SNESSetJacobian().
93362fef451SLois Curfman McInnes 
934c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
935c7afd0dbSLois Curfman McInnes 
93662fef451SLois Curfman McInnes    Input Parameters:
937c7afd0dbSLois Curfman McInnes +  snes - the SNES context
938c7afd0dbSLois Curfman McInnes -  x - input vector
93962fef451SLois Curfman McInnes 
94062fef451SLois Curfman McInnes    Output Parameters:
941c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
94262fef451SLois Curfman McInnes .  B - optional preconditioning matrix
943c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
944fee21e36SBarry Smith 
94562fef451SLois Curfman McInnes    Notes:
94662fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
94762fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
94862fef451SLois Curfman McInnes 
949dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
950dc5a77f8SLois Curfman McInnes    flag parameter.
95162fef451SLois Curfman McInnes 
95262fef451SLois Curfman McInnes    SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS
95362fef451SLois Curfman McInnes    methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION
954005c665bSBarry Smith    methods is SNESComputeHessian().
95562fef451SLois Curfman McInnes 
95636851e7fSLois Curfman McInnes    Level: developer
95736851e7fSLois Curfman McInnes 
95862fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
95962fef451SLois Curfman McInnes 
96062fef451SLois Curfman McInnes .seealso:  SNESSetJacobian(), SLESSetOperators()
96162fef451SLois Curfman McInnes @*/
9629b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
9639b94acceSBarry Smith {
9649b94acceSBarry Smith   int    ierr;
9653a40ed3dSBarry Smith 
9663a40ed3dSBarry Smith   PetscFunctionBegin;
967184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
968184914b5SBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE);
969184914b5SBarry Smith   PetscCheckSameComm(snes,X);
9703a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
971*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
9723a40ed3dSBarry Smith   }
9733a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
97481bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
975c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
976d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
9779b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
978d64ed03dSBarry Smith   PetscStackPop;
97981bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
9806d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
98177c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
98277c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
9833a40ed3dSBarry Smith   PetscFunctionReturn(0);
9849b94acceSBarry Smith }
9859b94acceSBarry Smith 
9865615d1e5SSatish Balay #undef __FUNC__
987b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESComputeHessian"
98862fef451SLois Curfman McInnes /*@
98962fef451SLois Curfman McInnes    SNESComputeHessian - Computes the Hessian matrix that has been
99062fef451SLois Curfman McInnes    set with SNESSetHessian().
99162fef451SLois Curfman McInnes 
992c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
993c7afd0dbSLois Curfman McInnes 
99462fef451SLois Curfman McInnes    Input Parameters:
995c7afd0dbSLois Curfman McInnes +  snes - the SNES context
996c7afd0dbSLois Curfman McInnes -  x - input vector
99762fef451SLois Curfman McInnes 
99862fef451SLois Curfman McInnes    Output Parameters:
999c7afd0dbSLois Curfman McInnes +  A - Hessian matrix
100062fef451SLois Curfman McInnes .  B - optional preconditioning matrix
1001c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
1002fee21e36SBarry Smith 
100362fef451SLois Curfman McInnes    Notes:
100462fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
100562fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
100662fef451SLois Curfman McInnes 
1007dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the
1008dc5a77f8SLois Curfman McInnes    flag parameter.
100962fef451SLois Curfman McInnes 
101062fef451SLois Curfman McInnes    SNESComputeHessian() is valid only for
101162fef451SLois Curfman McInnes    SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for
101262fef451SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian().
101362fef451SLois Curfman McInnes 
101436851e7fSLois Curfman McInnes    SNESComputeHessian() is typically used within minimization
101536851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
101636851e7fSLois Curfman McInnes    themselves.
101736851e7fSLois Curfman McInnes 
101836851e7fSLois Curfman McInnes    Level: developer
101936851e7fSLois Curfman McInnes 
102062fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix
102162fef451SLois Curfman McInnes 
1022a86d99e1SLois Curfman McInnes .seealso:  SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(),
1023a86d99e1SLois Curfman McInnes            SNESComputeMinimizationFunction()
102462fef451SLois Curfman McInnes @*/
102562fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag)
10269b94acceSBarry Smith {
10279b94acceSBarry Smith   int    ierr;
10283a40ed3dSBarry Smith 
10293a40ed3dSBarry Smith   PetscFunctionBegin;
1030184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1031184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
1032184914b5SBarry Smith   PetscCheckSameComm(snes,x);
10333a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1034*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
10353a40ed3dSBarry Smith   }
10363a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
103781bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr);
10380f4a323eSLois Curfman McInnes   *flag = DIFFERENT_NONZERO_PATTERN;
1039d64ed03dSBarry Smith   PetscStackPush("SNES user Hessian function");
104062fef451SLois Curfman McInnes   ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP);CHKERRQ(ierr);
1041d64ed03dSBarry Smith   PetscStackPop;
104281bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_HessianEval,snes,x,*A,*B);CHKERRQ(ierr);
10430f4a323eSLois Curfman McInnes   /* make sure user returned a correct Jacobian and preconditioner */
104477c4ece6SBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE);
104577c4ece6SBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE);
10463a40ed3dSBarry Smith   PetscFunctionReturn(0);
10479b94acceSBarry Smith }
10489b94acceSBarry Smith 
10495615d1e5SSatish Balay #undef __FUNC__
1050b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetJacobian"
10519b94acceSBarry Smith /*@C
10529b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
1053044dda88SLois Curfman McInnes    location to store the matrix.
10549b94acceSBarry Smith 
1055c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1056c7afd0dbSLois Curfman McInnes 
10579b94acceSBarry Smith    Input Parameters:
1058c7afd0dbSLois Curfman McInnes +  snes - the SNES context
10599b94acceSBarry Smith .  A - Jacobian matrix
10609b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
10619b94acceSBarry Smith .  func - Jacobian evaluation routine
1062c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
10632cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
10649b94acceSBarry Smith 
10659b94acceSBarry Smith    Calling sequence of func:
10668d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
10679b94acceSBarry Smith 
1068c7afd0dbSLois Curfman McInnes +  x - input vector
10699b94acceSBarry Smith .  A - Jacobian matrix
10709b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1071ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1072ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
1073c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
10749b94acceSBarry Smith 
10759b94acceSBarry Smith    Notes:
1076dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
10772cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1078ac21db08SLois Curfman McInnes 
1079ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
10809b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
10819b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
10829b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
10839b94acceSBarry Smith    throughout the global iterations.
10849b94acceSBarry Smith 
108536851e7fSLois Curfman McInnes    Level: beginner
108636851e7fSLois Curfman McInnes 
10879b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
10889b94acceSBarry Smith 
1089ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction()
10909b94acceSBarry Smith @*/
1091454a90a3SBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
10929b94acceSBarry Smith {
10933a40ed3dSBarry Smith   PetscFunctionBegin;
109477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
109500036973SBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE);
109600036973SBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE);
109700036973SBarry Smith   if (A) PetscCheckSameComm(snes,A);
109800036973SBarry Smith   if (B) PetscCheckSameComm(snes,B);
1099a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1100*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
1101a8c6a408SBarry Smith   }
1102184914b5SBarry Smith 
11039b94acceSBarry Smith   snes->computejacobian = func;
11049b94acceSBarry Smith   snes->jacP            = ctx;
11059b94acceSBarry Smith   snes->jacobian        = A;
11069b94acceSBarry Smith   snes->jacobian_pre    = B;
11073a40ed3dSBarry Smith   PetscFunctionReturn(0);
11089b94acceSBarry Smith }
110962fef451SLois Curfman McInnes 
11105615d1e5SSatish Balay #undef __FUNC__
1111b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetJacobian"
1112c2aafc4cSSatish Balay /*@C
1113b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
1114b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
1115b4fd4287SBarry Smith 
1116c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
1117c7afd0dbSLois Curfman McInnes 
1118b4fd4287SBarry Smith    Input Parameter:
1119b4fd4287SBarry Smith .  snes - the nonlinear solver context
1120b4fd4287SBarry Smith 
1121b4fd4287SBarry Smith    Output Parameters:
1122c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
1123b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
112400036973SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
112500036973SBarry Smith -  func - location to put Jacobian function (or PETSC_NULL)
1126fee21e36SBarry Smith 
112736851e7fSLois Curfman McInnes    Level: advanced
112836851e7fSLois Curfman McInnes 
1129b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
1130b4fd4287SBarry Smith @*/
113100036973SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
1132b4fd4287SBarry Smith {
11333a40ed3dSBarry Smith   PetscFunctionBegin;
1134184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
11353a40ed3dSBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
1136*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
11373a40ed3dSBarry Smith   }
1138b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
1139b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
1140b4fd4287SBarry Smith   if (ctx)  *ctx  = snes->jacP;
114100036973SBarry Smith   if (func) *func = snes->computejacobian;
11423a40ed3dSBarry Smith   PetscFunctionReturn(0);
1143b4fd4287SBarry Smith }
1144b4fd4287SBarry Smith 
11455615d1e5SSatish Balay #undef __FUNC__
1146b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetHessian"
11479b94acceSBarry Smith /*@C
11489b94acceSBarry Smith    SNESSetHessian - Sets the function to compute Hessian as well as the
1149044dda88SLois Curfman McInnes    location to store the matrix.
11509b94acceSBarry Smith 
1151c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
1152c7afd0dbSLois Curfman McInnes 
11539b94acceSBarry Smith    Input Parameters:
1154c7afd0dbSLois Curfman McInnes +  snes - the SNES context
11559b94acceSBarry Smith .  A - Hessian matrix
11569b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Hessian)
11579b94acceSBarry Smith .  func - Jacobian evaluation routine
1158c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
1159313e4042SLois Curfman McInnes          Hessian evaluation routine (may be PETSC_NULL)
11609b94acceSBarry Smith 
11619b94acceSBarry Smith    Calling sequence of func:
11628d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
11639b94acceSBarry Smith 
1164c7afd0dbSLois Curfman McInnes +  x - input vector
11659b94acceSBarry Smith .  A - Hessian matrix
11669b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
1167ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
1168ac21db08SLois Curfman McInnes    structure (same as flag in SLESSetOperators())
1169c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Hessian context
11709b94acceSBarry Smith 
11719b94acceSBarry Smith    Notes:
1172dc5a77f8SLois Curfman McInnes    See SLESSetOperators() for important information about setting the flag
11732cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
1174ac21db08SLois Curfman McInnes 
11759b94acceSBarry Smith    The function func() takes Mat * as the matrix arguments rather than Mat.
11769b94acceSBarry Smith    This allows the Hessian evaluation routine to replace A and/or B with a
11779b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
11789b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
11799b94acceSBarry Smith    throughout the global iterations.
11809b94acceSBarry Smith 
118136851e7fSLois Curfman McInnes    Level: beginner
118236851e7fSLois Curfman McInnes 
11839b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix
11849b94acceSBarry Smith 
1185ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators()
11869b94acceSBarry Smith @*/
1187454a90a3SBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
11889b94acceSBarry Smith {
11893a40ed3dSBarry Smith   PetscFunctionBegin;
119077c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1191184914b5SBarry Smith   PetscValidHeaderSpecific(A,MAT_COOKIE);
1192184914b5SBarry Smith   PetscValidHeaderSpecific(B,MAT_COOKIE);
1193184914b5SBarry Smith   PetscCheckSameComm(snes,A);
1194184914b5SBarry Smith   PetscCheckSameComm(snes,B);
1195d4bb536fSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
1196*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
1197d4bb536fSBarry Smith   }
11989b94acceSBarry Smith   snes->computejacobian = func;
11999b94acceSBarry Smith   snes->jacP            = ctx;
12009b94acceSBarry Smith   snes->jacobian        = A;
12019b94acceSBarry Smith   snes->jacobian_pre    = B;
12023a40ed3dSBarry Smith   PetscFunctionReturn(0);
12039b94acceSBarry Smith }
12049b94acceSBarry Smith 
12055615d1e5SSatish Balay #undef __FUNC__
1206b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetHessian"
120762fef451SLois Curfman McInnes /*@
120862fef451SLois Curfman McInnes    SNESGetHessian - Returns the Hessian matrix and optionally the user
120962fef451SLois Curfman McInnes    provided context for evaluating the Hessian.
121062fef451SLois Curfman McInnes 
1211c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object is parallel if SNES object is parallel
1212c7afd0dbSLois Curfman McInnes 
121362fef451SLois Curfman McInnes    Input Parameter:
121462fef451SLois Curfman McInnes .  snes - the nonlinear solver context
121562fef451SLois Curfman McInnes 
121662fef451SLois Curfman McInnes    Output Parameters:
1217c7afd0dbSLois Curfman McInnes +  A - location to stash Hessian matrix (or PETSC_NULL)
121862fef451SLois Curfman McInnes .  B - location to stash preconditioner matrix (or PETSC_NULL)
1219c7afd0dbSLois Curfman McInnes -  ctx - location to stash Hessian ctx (or PETSC_NULL)
1220fee21e36SBarry Smith 
122136851e7fSLois Curfman McInnes    Level: advanced
122236851e7fSLois Curfman McInnes 
122362fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian()
1224c7afd0dbSLois Curfman McInnes 
1225c7afd0dbSLois Curfman McInnes .keywords: SNES, get, Hessian
122662fef451SLois Curfman McInnes @*/
122762fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B,void **ctx)
122862fef451SLois Curfman McInnes {
12293a40ed3dSBarry Smith   PetscFunctionBegin;
1230184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
12313a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){
1232*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
12333a40ed3dSBarry Smith   }
123462fef451SLois Curfman McInnes   if (A)   *A = snes->jacobian;
123562fef451SLois Curfman McInnes   if (B)   *B = snes->jacobian_pre;
123662fef451SLois Curfman McInnes   if (ctx) *ctx = snes->jacP;
12373a40ed3dSBarry Smith   PetscFunctionReturn(0);
123862fef451SLois Curfman McInnes }
123962fef451SLois Curfman McInnes 
12409b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
12419b94acceSBarry Smith 
12425615d1e5SSatish Balay #undef __FUNC__
1243b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetUp"
12449b94acceSBarry Smith /*@
12459b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
1246272ac6f2SLois Curfman McInnes    of a nonlinear solver.
12479b94acceSBarry Smith 
1248fee21e36SBarry Smith    Collective on SNES
1249fee21e36SBarry Smith 
1250c7afd0dbSLois Curfman McInnes    Input Parameters:
1251c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1252c7afd0dbSLois Curfman McInnes -  x - the solution vector
1253c7afd0dbSLois Curfman McInnes 
1254272ac6f2SLois Curfman McInnes    Notes:
1255272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
1256272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
1257272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
1258272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
1259272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
1260272ac6f2SLois Curfman McInnes 
126136851e7fSLois Curfman McInnes    Level: advanced
126236851e7fSLois Curfman McInnes 
12639b94acceSBarry Smith .keywords: SNES, nonlinear, setup
12649b94acceSBarry Smith 
12659b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
12669b94acceSBarry Smith @*/
12678ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
12689b94acceSBarry Smith {
1269f1af5d2fSBarry Smith   int        ierr;
1270f1af5d2fSBarry Smith   PetscTruth flg;
12713a40ed3dSBarry Smith 
12723a40ed3dSBarry Smith   PetscFunctionBegin;
127377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
127477c4ece6SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
1275184914b5SBarry Smith   PetscCheckSameComm(snes,x);
12768ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
12779b94acceSBarry Smith 
1278c1f60f51SBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
1279c1f60f51SBarry Smith   /*
1280c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
1281dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
1282c1f60f51SBarry Smith   */
1283c1f60f51SBarry Smith   if (flg) {
1284c1f60f51SBarry Smith     Mat J;
12855a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1286c1f60f51SBarry Smith     PLogObjectParent(snes,J);
1287c1f60f51SBarry Smith     snes->mfshell = J;
1288c1f60f51SBarry Smith     snes->jacobian = J;
1289c2aafc4cSSatish Balay     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
1290c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n");
1291d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
1292c1f60f51SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n");
1293d4bb536fSBarry Smith     } else {
1294*29bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Method class doesn't support matrix-free operator option");
1295d4bb536fSBarry Smith     }
12965a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1297c1f60f51SBarry Smith   }
1298272ac6f2SLois Curfman McInnes   ierr = OptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1299c1f60f51SBarry Smith   /*
1300dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1301c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1302c1f60f51SBarry Smith    */
130331615d04SLois Curfman McInnes   if (flg) {
1304272ac6f2SLois Curfman McInnes     Mat J;
13055a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
1306272ac6f2SLois Curfman McInnes     PLogObjectParent(snes,J);
1307272ac6f2SLois Curfman McInnes     snes->mfshell = J;
130883e56afdSLois Curfman McInnes     if (snes->method_class == SNES_NONLINEAR_EQUATIONS) {
13091d1367b7SBarry Smith       ierr = SNESSetJacobian(snes,J,J,MatSNESMFFormJacobian,snes->funP);CHKERRQ(ierr);
131094a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n");
1311d64ed03dSBarry Smith     } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) {
13121d1367b7SBarry Smith       ierr = SNESSetHessian(snes,J,J,MatSNESMFFormJacobian,snes->funP);CHKERRQ(ierr);
131394a424c1SBarry Smith       PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n");
1314d4bb536fSBarry Smith     } else {
1315*29bbc08cSBarry Smith       SETERRQ(PETSC_ERR_SUP,"Method class doesn't support matrix-free option");
1316d4bb536fSBarry Smith     }
13175a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
1318272ac6f2SLois Curfman McInnes   }
13199b94acceSBarry Smith   if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) {
13206831982aSBarry Smith     PetscTruth iseqtr;
13216831982aSBarry Smith 
1322*29bbc08cSBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1323*29bbc08cSBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
1324*29bbc08cSBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1325a8c6a408SBarry Smith     if (snes->vec_func == snes->vec_sol) {
1326*29bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1327a8c6a408SBarry Smith     }
1328a703fe33SLois Curfman McInnes 
1329a703fe33SLois Curfman McInnes     /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
13306831982aSBarry Smith     ierr = PetscTypeCompare((PetscObject)snes,SNESEQTR,&iseqtr);CHKERRQ(ierr);
13316831982aSBarry Smith     if (snes->ksp_ewconv && !iseqtr) {
1332a703fe33SLois Curfman McInnes       SLES sles; KSP ksp;
1333a703fe33SLois Curfman McInnes       ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr);
1334a703fe33SLois Curfman McInnes       ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr);
1335186905e3SBarry Smith       ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1336a703fe33SLois Curfman McInnes     }
13379b94acceSBarry Smith   } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) {
1338*29bbc08cSBarry Smith     if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first");
1339*29bbc08cSBarry Smith     if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetGradient() first");
1340a8c6a408SBarry Smith     if (!snes->computeumfunction) {
1341*29bbc08cSBarry Smith       SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetMinimizationFunction() first");
1342a8c6a408SBarry Smith     }
1343*29bbc08cSBarry Smith     if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetHessian()");
1344d4bb536fSBarry Smith   } else {
1345*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Unknown method class");
1346d4bb536fSBarry Smith   }
1347a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
134882bf6240SBarry Smith   snes->setupcalled = 1;
13493a40ed3dSBarry Smith   PetscFunctionReturn(0);
13509b94acceSBarry Smith }
13519b94acceSBarry Smith 
13525615d1e5SSatish Balay #undef __FUNC__
1353b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESDestroy"
13549b94acceSBarry Smith /*@C
13559b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
13569b94acceSBarry Smith    with SNESCreate().
13579b94acceSBarry Smith 
1358c7afd0dbSLois Curfman McInnes    Collective on SNES
1359c7afd0dbSLois Curfman McInnes 
13609b94acceSBarry Smith    Input Parameter:
13619b94acceSBarry Smith .  snes - the SNES context
13629b94acceSBarry Smith 
136336851e7fSLois Curfman McInnes    Level: beginner
136436851e7fSLois Curfman McInnes 
13659b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
13669b94acceSBarry Smith 
136763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
13689b94acceSBarry Smith @*/
13699b94acceSBarry Smith int SNESDestroy(SNES snes)
13709b94acceSBarry Smith {
1371b8a78c4aSBarry Smith   int i,ierr;
13723a40ed3dSBarry Smith 
13733a40ed3dSBarry Smith   PetscFunctionBegin;
137477c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
13753a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1376d4bb536fSBarry Smith 
1377be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
13780f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1379be0abb6dSBarry Smith 
1380e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1381606d414cSSatish Balay   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
1382522c5e43SBarry Smith   if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);}
13839b94acceSBarry Smith   ierr = SLESDestroy(snes->sles);CHKERRQ(ierr);
1384522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1385b8a78c4aSBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1386b8a78c4aSBarry Smith     if (snes->monitordestroy[i]) {
1387b8a78c4aSBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1388b8a78c4aSBarry Smith     }
1389b8a78c4aSBarry Smith   }
13909b94acceSBarry Smith   PLogObjectDestroy((PetscObject)snes);
13910452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
13923a40ed3dSBarry Smith   PetscFunctionReturn(0);
13939b94acceSBarry Smith }
13949b94acceSBarry Smith 
13959b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
13969b94acceSBarry Smith 
13975615d1e5SSatish Balay #undef __FUNC__
1398b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetTolerances"
13999b94acceSBarry Smith /*@
1400d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
14019b94acceSBarry Smith 
1402c7afd0dbSLois Curfman McInnes    Collective on SNES
1403c7afd0dbSLois Curfman McInnes 
14049b94acceSBarry Smith    Input Parameters:
1405c7afd0dbSLois Curfman McInnes +  snes - the SNES context
140633174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
140733174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
140833174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
140933174efeSLois Curfman McInnes            of the change in the solution between steps
141033174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1411c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1412fee21e36SBarry Smith 
141333174efeSLois Curfman McInnes    Options Database Keys:
1414c7afd0dbSLois Curfman McInnes +    -snes_atol <atol> - Sets atol
1415c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1416c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1417c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1418c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
14199b94acceSBarry Smith 
1420d7a720efSLois Curfman McInnes    Notes:
14219b94acceSBarry Smith    The default maximum number of iterations is 50.
14229b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
14239b94acceSBarry Smith 
142436851e7fSLois Curfman McInnes    Level: intermediate
142536851e7fSLois Curfman McInnes 
142633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
14279b94acceSBarry Smith 
1428d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
14299b94acceSBarry Smith @*/
1430329f5518SBarry Smith int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
14319b94acceSBarry Smith {
14323a40ed3dSBarry Smith   PetscFunctionBegin;
143377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1434d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1435d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1436d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1437d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1438d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
14393a40ed3dSBarry Smith   PetscFunctionReturn(0);
14409b94acceSBarry Smith }
14419b94acceSBarry Smith 
14425615d1e5SSatish Balay #undef __FUNC__
1443b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetTolerances"
14449b94acceSBarry Smith /*@
144533174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
144633174efeSLois Curfman McInnes 
1447c7afd0dbSLois Curfman McInnes    Not Collective
1448c7afd0dbSLois Curfman McInnes 
144933174efeSLois Curfman McInnes    Input Parameters:
1450c7afd0dbSLois Curfman McInnes +  snes - the SNES context
145133174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
145233174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
145333174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
145433174efeSLois Curfman McInnes            of the change in the solution between steps
145533174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1456c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1457fee21e36SBarry Smith 
145833174efeSLois Curfman McInnes    Notes:
145933174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
146033174efeSLois Curfman McInnes 
146136851e7fSLois Curfman McInnes    Level: intermediate
146236851e7fSLois Curfman McInnes 
146333174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
146433174efeSLois Curfman McInnes 
146533174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
146633174efeSLois Curfman McInnes @*/
1467329f5518SBarry Smith int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
146833174efeSLois Curfman McInnes {
14693a40ed3dSBarry Smith   PetscFunctionBegin;
147033174efeSLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
147133174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
147233174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
147333174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
147433174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
147533174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
14763a40ed3dSBarry Smith   PetscFunctionReturn(0);
147733174efeSLois Curfman McInnes }
147833174efeSLois Curfman McInnes 
14795615d1e5SSatish Balay #undef __FUNC__
1480b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetTrustRegionTolerance"
148133174efeSLois Curfman McInnes /*@
14829b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
14839b94acceSBarry Smith 
1484fee21e36SBarry Smith    Collective on SNES
1485fee21e36SBarry Smith 
1486c7afd0dbSLois Curfman McInnes    Input Parameters:
1487c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1488c7afd0dbSLois Curfman McInnes -  tol - tolerance
1489c7afd0dbSLois Curfman McInnes 
14909b94acceSBarry Smith    Options Database Key:
1491c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
14929b94acceSBarry Smith 
149336851e7fSLois Curfman McInnes    Level: intermediate
149436851e7fSLois Curfman McInnes 
14959b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
14969b94acceSBarry Smith 
1497d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
14989b94acceSBarry Smith @*/
1499329f5518SBarry Smith int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
15009b94acceSBarry Smith {
15013a40ed3dSBarry Smith   PetscFunctionBegin;
150277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15039b94acceSBarry Smith   snes->deltatol = tol;
15043a40ed3dSBarry Smith   PetscFunctionReturn(0);
15059b94acceSBarry Smith }
15069b94acceSBarry Smith 
15075615d1e5SSatish Balay #undef __FUNC__
1508b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetMinimizationFunctionTolerance"
15099b94acceSBarry Smith /*@
151077c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance
15119b94acceSBarry Smith    for unconstrained minimization solvers.
15129b94acceSBarry Smith 
1513fee21e36SBarry Smith    Collective on SNES
1514fee21e36SBarry Smith 
1515c7afd0dbSLois Curfman McInnes    Input Parameters:
1516c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1517c7afd0dbSLois Curfman McInnes -  ftol - minimum function tolerance
1518c7afd0dbSLois Curfman McInnes 
15199b94acceSBarry Smith    Options Database Key:
1520c7afd0dbSLois Curfman McInnes .  -snes_fmin <ftol> - Sets ftol
15219b94acceSBarry Smith 
15229b94acceSBarry Smith    Note:
152377c4ece6SBarry Smith    SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION
15249b94acceSBarry Smith    methods only.
15259b94acceSBarry Smith 
152636851e7fSLois Curfman McInnes    Level: intermediate
152736851e7fSLois Curfman McInnes 
15289b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance
15299b94acceSBarry Smith 
1530d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance()
15319b94acceSBarry Smith @*/
1532329f5518SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,PetscReal ftol)
15339b94acceSBarry Smith {
15343a40ed3dSBarry Smith   PetscFunctionBegin;
153577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
15369b94acceSBarry Smith   snes->fmin = ftol;
15373a40ed3dSBarry Smith   PetscFunctionReturn(0);
15389b94acceSBarry Smith }
1539df9fa365SBarry Smith /*
1540df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1541df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1542df9fa365SBarry Smith    macros instead of functions
1543df9fa365SBarry Smith */
1544ce1608b8SBarry Smith #undef __FUNC__
1545b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESLGMonitor"
1546329f5518SBarry Smith int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1547ce1608b8SBarry Smith {
1548ce1608b8SBarry Smith   int ierr;
1549ce1608b8SBarry Smith 
1550ce1608b8SBarry Smith   PetscFunctionBegin;
1551184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1552ce1608b8SBarry Smith   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1553ce1608b8SBarry Smith   PetscFunctionReturn(0);
1554ce1608b8SBarry Smith }
1555ce1608b8SBarry Smith 
1556df9fa365SBarry Smith #undef __FUNC__
1557b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESLGMonitorCreate"
1558df9fa365SBarry Smith int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,DrawLG *draw)
1559df9fa365SBarry Smith {
1560df9fa365SBarry Smith   int ierr;
1561df9fa365SBarry Smith 
1562df9fa365SBarry Smith   PetscFunctionBegin;
1563df9fa365SBarry Smith   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1564df9fa365SBarry Smith   PetscFunctionReturn(0);
1565df9fa365SBarry Smith }
1566df9fa365SBarry Smith 
1567df9fa365SBarry Smith #undef __FUNC__
1568b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESLGMonitorDestroy"
1569df9fa365SBarry Smith int SNESLGMonitorDestroy(DrawLG draw)
1570df9fa365SBarry Smith {
1571df9fa365SBarry Smith   int ierr;
1572df9fa365SBarry Smith 
1573df9fa365SBarry Smith   PetscFunctionBegin;
1574df9fa365SBarry Smith   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1575df9fa365SBarry Smith   PetscFunctionReturn(0);
1576df9fa365SBarry Smith }
1577df9fa365SBarry Smith 
15789b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
15799b94acceSBarry Smith 
15805615d1e5SSatish Balay #undef __FUNC__
1581b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetMonitor"
15829b94acceSBarry Smith /*@C
15835cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
15849b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
15859b94acceSBarry Smith    progress.
15869b94acceSBarry Smith 
1587fee21e36SBarry Smith    Collective on SNES
1588fee21e36SBarry Smith 
1589c7afd0dbSLois Curfman McInnes    Input Parameters:
1590c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1591c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1592b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1593b3006f0bSLois Curfman McInnes           monitor routine (use PETSC_NULL if no context is desitre)
1594b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1595b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
15969b94acceSBarry Smith 
1597c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1598329f5518SBarry Smith $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1599c7afd0dbSLois Curfman McInnes 
1600c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1601c7afd0dbSLois Curfman McInnes .    its - iteration number
1602c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
1603c7afd0dbSLois Curfman McInnes             for SNES_NONLINEAR_EQUATIONS methods
160440a0c1c6SLois Curfman McInnes .    norm - 2-norm gradient value (may be estimated)
1605c7afd0dbSLois Curfman McInnes             for SNES_UNCONSTRAINED_MINIMIZATION methods
160640a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
16079b94acceSBarry Smith 
16089665c990SLois Curfman McInnes    Options Database Keys:
1609c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1610c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1611c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1612c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1613c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1614c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1615c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1616c7afd0dbSLois Curfman McInnes                             the options database.
16179665c990SLois Curfman McInnes 
1618639f9d9dSBarry Smith    Notes:
16196bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
16206bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
16216bc08f3fSLois Curfman McInnes    order in which they were set.
1622639f9d9dSBarry Smith 
162336851e7fSLois Curfman McInnes    Level: intermediate
162436851e7fSLois Curfman McInnes 
16259b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
16269b94acceSBarry Smith 
16275cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
16289b94acceSBarry Smith @*/
1629329f5518SBarry Smith int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *))
16309b94acceSBarry Smith {
16313a40ed3dSBarry Smith   PetscFunctionBegin;
1632184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1633639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
1634*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1635639f9d9dSBarry Smith   }
1636639f9d9dSBarry Smith 
1637639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1638b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1639639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
16403a40ed3dSBarry Smith   PetscFunctionReturn(0);
16419b94acceSBarry Smith }
16429b94acceSBarry Smith 
16435615d1e5SSatish Balay #undef __FUNC__
1644b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESClearMonitor"
16455cd90555SBarry Smith /*@C
16465cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
16475cd90555SBarry Smith 
1648c7afd0dbSLois Curfman McInnes    Collective on SNES
1649c7afd0dbSLois Curfman McInnes 
16505cd90555SBarry Smith    Input Parameters:
16515cd90555SBarry Smith .  snes - the SNES context
16525cd90555SBarry Smith 
16535cd90555SBarry Smith    Options Database:
1654c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1655c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1656c7afd0dbSLois Curfman McInnes     set via the options database
16575cd90555SBarry Smith 
16585cd90555SBarry Smith    Notes:
16595cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
16605cd90555SBarry Smith 
166136851e7fSLois Curfman McInnes    Level: intermediate
166236851e7fSLois Curfman McInnes 
16635cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
16645cd90555SBarry Smith 
16655cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
16665cd90555SBarry Smith @*/
16675cd90555SBarry Smith int SNESClearMonitor(SNES snes)
16685cd90555SBarry Smith {
16695cd90555SBarry Smith   PetscFunctionBegin;
1670184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
16715cd90555SBarry Smith   snes->numbermonitors = 0;
16725cd90555SBarry Smith   PetscFunctionReturn(0);
16735cd90555SBarry Smith }
16745cd90555SBarry Smith 
16755cd90555SBarry Smith #undef __FUNC__
1676b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetConvergenceTest"
16779b94acceSBarry Smith /*@C
16789b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
16799b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
16809b94acceSBarry Smith 
1681fee21e36SBarry Smith    Collective on SNES
1682fee21e36SBarry Smith 
1683c7afd0dbSLois Curfman McInnes    Input Parameters:
1684c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1685c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1686c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1687c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
16889b94acceSBarry Smith 
1689c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1690329f5518SBarry Smith $     int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1691c7afd0dbSLois Curfman McInnes 
1692c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1693c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1694184914b5SBarry Smith .    reason - reason for convergence/divergence
1695c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
1696c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods)
1697c7afd0dbSLois Curfman McInnes .    f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods)
1698c7afd0dbSLois Curfman McInnes .    gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods)
1699c7afd0dbSLois Curfman McInnes -    f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods)
17009b94acceSBarry Smith 
170136851e7fSLois Curfman McInnes    Level: advanced
170236851e7fSLois Curfman McInnes 
17039b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
17049b94acceSBarry Smith 
170540191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
170640191667SLois Curfman McInnes           SNESConverged_UM_LS(), SNESConverged_UM_TR()
17079b94acceSBarry Smith @*/
1708329f5518SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
17099b94acceSBarry Smith {
17103a40ed3dSBarry Smith   PetscFunctionBegin;
1711184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
17129b94acceSBarry Smith   (snes)->converged = func;
17139b94acceSBarry Smith   (snes)->cnvP      = cctx;
17143a40ed3dSBarry Smith   PetscFunctionReturn(0);
17159b94acceSBarry Smith }
17169b94acceSBarry Smith 
17175615d1e5SSatish Balay #undef __FUNC__
1718b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetConvergedReason"
1719184914b5SBarry Smith /*@C
1720184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1721184914b5SBarry Smith 
1722184914b5SBarry Smith    Not Collective
1723184914b5SBarry Smith 
1724184914b5SBarry Smith    Input Parameter:
1725184914b5SBarry Smith .  snes - the SNES context
1726184914b5SBarry Smith 
1727184914b5SBarry Smith    Output Parameter:
1728e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1729184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1730184914b5SBarry Smith 
1731184914b5SBarry Smith    Level: intermediate
1732184914b5SBarry Smith 
1733184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1734184914b5SBarry Smith 
1735184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1736184914b5SBarry Smith 
1737184914b5SBarry Smith .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(),
1738184914b5SBarry Smith           SNESConverged_UM_LS(), SNESConverged_UM_TR()
1739184914b5SBarry Smith @*/
1740184914b5SBarry Smith int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1741184914b5SBarry Smith {
1742184914b5SBarry Smith   PetscFunctionBegin;
1743184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1744184914b5SBarry Smith   *reason = snes->reason;
1745184914b5SBarry Smith   PetscFunctionReturn(0);
1746184914b5SBarry Smith }
1747184914b5SBarry Smith 
1748184914b5SBarry Smith #undef __FUNC__
1749b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetConvergenceHistory"
1750c9005455SLois Curfman McInnes /*@
1751c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1752c9005455SLois Curfman McInnes 
1753fee21e36SBarry Smith    Collective on SNES
1754fee21e36SBarry Smith 
1755c7afd0dbSLois Curfman McInnes    Input Parameters:
1756c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1757c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1758cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1759758f92a0SBarry Smith .  na  - size of a and its
1760758f92a0SBarry Smith -  reset - PETSC_TRUTH indicates each new nonlinear solve resets the history counter to zero,
1761758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1762c7afd0dbSLois Curfman McInnes 
1763c9005455SLois Curfman McInnes    Notes:
1764c9005455SLois Curfman McInnes    If set, this array will contain the function norms (for
1765c9005455SLois Curfman McInnes    SNES_NONLINEAR_EQUATIONS methods) or gradient norms
1766c9005455SLois Curfman McInnes    (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed
1767c9005455SLois Curfman McInnes    at each step.
1768c9005455SLois Curfman McInnes 
1769c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1770c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1771c9005455SLois Curfman McInnes    during the section of code that is being timed.
1772c9005455SLois Curfman McInnes 
177336851e7fSLois Curfman McInnes    Level: intermediate
177436851e7fSLois Curfman McInnes 
1775c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1776758f92a0SBarry Smith 
177708405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1778758f92a0SBarry Smith 
1779c9005455SLois Curfman McInnes @*/
1780329f5518SBarry Smith int SNESSetConvergenceHistory(SNES snes,PetscReal *a,int *its,int na,PetscTruth reset)
1781c9005455SLois Curfman McInnes {
17823a40ed3dSBarry Smith   PetscFunctionBegin;
1783c9005455SLois Curfman McInnes   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1784c9005455SLois Curfman McInnes   if (na) PetscValidScalarPointer(a);
1785c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1786758f92a0SBarry Smith   snes->conv_hist_its   = its;
1787758f92a0SBarry Smith   snes->conv_hist_max   = na;
1788758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1789758f92a0SBarry Smith   PetscFunctionReturn(0);
1790758f92a0SBarry Smith }
1791758f92a0SBarry Smith 
1792758f92a0SBarry Smith #undef __FUNC__
1793b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetConvergenceHistory"
17940c4c9dddSBarry Smith /*@C
1795758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1796758f92a0SBarry Smith 
1797758f92a0SBarry Smith    Collective on SNES
1798758f92a0SBarry Smith 
1799758f92a0SBarry Smith    Input Parameter:
1800758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1801758f92a0SBarry Smith 
1802758f92a0SBarry Smith    Output Parameters:
1803758f92a0SBarry Smith .  a   - array to hold history
1804758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1805758f92a0SBarry Smith          negative if not converged) for each solve.
1806758f92a0SBarry Smith -  na  - size of a and its
1807758f92a0SBarry Smith 
1808758f92a0SBarry Smith    Notes:
1809758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1810758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1811758f92a0SBarry Smith 
1812758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1813758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1814758f92a0SBarry Smith    during the section of code that is being timed.
1815758f92a0SBarry Smith 
1816758f92a0SBarry Smith    Level: intermediate
1817758f92a0SBarry Smith 
1818758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1819758f92a0SBarry Smith 
1820758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1821758f92a0SBarry Smith 
1822758f92a0SBarry Smith @*/
1823329f5518SBarry Smith int SNESGetConvergenceHistory(SNES snes,PetscReal **a,int **its,int *na)
1824758f92a0SBarry Smith {
1825758f92a0SBarry Smith   PetscFunctionBegin;
1826758f92a0SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1827758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1828758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1829758f92a0SBarry Smith   if (na) *na   = snes->conv_hist_len;
18303a40ed3dSBarry Smith   PetscFunctionReturn(0);
1831c9005455SLois Curfman McInnes }
1832c9005455SLois Curfman McInnes 
1833c9005455SLois Curfman McInnes #undef __FUNC__
1834b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESScaleStep_Private"
18359b94acceSBarry Smith /*
18369b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
18379b94acceSBarry Smith    positive parameter delta.
18389b94acceSBarry Smith 
18399b94acceSBarry Smith     Input Parameters:
1840c7afd0dbSLois Curfman McInnes +   snes - the SNES context
18419b94acceSBarry Smith .   y - approximate solution of linear system
18429b94acceSBarry Smith .   fnorm - 2-norm of current function
1843c7afd0dbSLois Curfman McInnes -   delta - trust region size
18449b94acceSBarry Smith 
18459b94acceSBarry Smith     Output Parameters:
1846c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
18479b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
18489b94acceSBarry Smith     region, and exceeds zero otherwise.
1849c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
18509b94acceSBarry Smith 
18519b94acceSBarry Smith     Note:
18526831982aSBarry Smith     For non-trust region methods such as SNESEQLS, the parameter delta
18539b94acceSBarry Smith     is set to be the maximum allowable step size.
18549b94acceSBarry Smith 
18559b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
18569b94acceSBarry Smith */
1857329f5518SBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,
1858329f5518SBarry Smith                           PetscReal *gpnorm,PetscReal *ynorm)
18599b94acceSBarry Smith {
1860329f5518SBarry Smith   PetscReal norm;
18619b94acceSBarry Smith   Scalar cnorm;
18623a40ed3dSBarry Smith   int    ierr;
18633a40ed3dSBarry Smith 
18643a40ed3dSBarry Smith   PetscFunctionBegin;
1865184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1866184914b5SBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE);
1867184914b5SBarry Smith   PetscCheckSameComm(snes,y);
1868184914b5SBarry Smith 
18693a40ed3dSBarry Smith   ierr = VecNorm(y,NORM_2,&norm);CHKERRQ(ierr);
18709b94acceSBarry Smith   if (norm > *delta) {
18719b94acceSBarry Smith      norm = *delta/norm;
18729b94acceSBarry Smith      *gpnorm = (1.0 - norm)*(*fnorm);
18739b94acceSBarry Smith      cnorm = norm;
1874329f5518SBarry Smith      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
18759b94acceSBarry Smith      *ynorm = *delta;
18769b94acceSBarry Smith   } else {
18779b94acceSBarry Smith      *gpnorm = 0.0;
18789b94acceSBarry Smith      *ynorm = norm;
18799b94acceSBarry Smith   }
18803a40ed3dSBarry Smith   PetscFunctionReturn(0);
18819b94acceSBarry Smith }
18829b94acceSBarry Smith 
18835615d1e5SSatish Balay #undef __FUNC__
1884b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSolve"
18859b94acceSBarry Smith /*@
18869b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
188763a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
18889b94acceSBarry Smith 
1889c7afd0dbSLois Curfman McInnes    Collective on SNES
1890c7afd0dbSLois Curfman McInnes 
1891b2002411SLois Curfman McInnes    Input Parameters:
1892c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1893c7afd0dbSLois Curfman McInnes -  x - the solution vector
18949b94acceSBarry Smith 
18959b94acceSBarry Smith    Output Parameter:
1896b2002411SLois Curfman McInnes .  its - number of iterations until termination
18979b94acceSBarry Smith 
1898b2002411SLois Curfman McInnes    Notes:
18998ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
19008ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
19018ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
19028ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
19038ddd3da0SLois Curfman McInnes 
190436851e7fSLois Curfman McInnes    Level: beginner
190536851e7fSLois Curfman McInnes 
19069b94acceSBarry Smith .keywords: SNES, nonlinear, solve
19079b94acceSBarry Smith 
190863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy()
19099b94acceSBarry Smith @*/
19108ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its)
19119b94acceSBarry Smith {
1912f1af5d2fSBarry Smith   int        ierr;
1913f1af5d2fSBarry Smith   PetscTruth flg;
1914052efed2SBarry Smith 
19153a40ed3dSBarry Smith   PetscFunctionBegin;
191677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
1917184914b5SBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE);
1918184914b5SBarry Smith   PetscCheckSameComm(snes,x);
191974679c65SBarry Smith   PetscValidIntPointer(its);
1920*29bbc08cSBarry Smith   if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
192174637425SBarry Smith 
192282bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
1923c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
1924758f92a0SBarry Smith   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
192581bfdfe8SSatish Balay   ierr = PLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1926c96a6f78SLois Curfman McInnes   snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0;
19279b94acceSBarry Smith   ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr);
192881bfdfe8SSatish Balay   ierr = PLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
19290f5bd95cSBarry Smith   ierr = OptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
19306d4a8577SBarry Smith   if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD);CHKERRQ(ierr); }
19313a40ed3dSBarry Smith   PetscFunctionReturn(0);
19329b94acceSBarry Smith }
19339b94acceSBarry Smith 
19349b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
19359b94acceSBarry Smith 
19365615d1e5SSatish Balay #undef __FUNC__
1937b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetType"
193882bf6240SBarry Smith /*@C
19394b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
19409b94acceSBarry Smith 
1941fee21e36SBarry Smith    Collective on SNES
1942fee21e36SBarry Smith 
1943c7afd0dbSLois Curfman McInnes    Input Parameters:
1944c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1945454a90a3SBarry Smith -  type - a known method
1946c7afd0dbSLois Curfman McInnes 
1947c7afd0dbSLois Curfman McInnes    Options Database Key:
1948454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
1949c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1950ae12b187SLois Curfman McInnes 
19519b94acceSBarry Smith    Notes:
1952e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
19536831982aSBarry Smith +    SNESEQLS - Newton's method with line search
1954c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
19556831982aSBarry Smith .    SNESEQTR - Newton's method with trust region
1956c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
19576831982aSBarry Smith .    SNESUMTR - Newton's method with trust region
1958c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
19596831982aSBarry Smith -    SNESUMLS - Newton's method with line search
1960c7afd0dbSLois Curfman McInnes      (unconstrained minimization)
19619b94acceSBarry Smith 
1962ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1963ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1964ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1965ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1966ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1967ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1968ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1969ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1970ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
197136851e7fSLois Curfman McInnes   appropriate method.  In other words, this routine is not for beginners.
197236851e7fSLois Curfman McInnes 
197336851e7fSLois Curfman McInnes   Level: intermediate
1974a703fe33SLois Curfman McInnes 
1975454a90a3SBarry Smith .keywords: SNES, set, type
19769b94acceSBarry Smith @*/
1977454a90a3SBarry Smith int SNESSetType(SNES snes,SNESType type)
19789b94acceSBarry Smith {
19790f5bd95cSBarry Smith   int        ierr,(*r)(SNES);
19806831982aSBarry Smith   PetscTruth match;
19813a40ed3dSBarry Smith 
19823a40ed3dSBarry Smith   PetscFunctionBegin;
198377c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
19840f5bd95cSBarry Smith   PetscValidCharPointer(type);
198582bf6240SBarry Smith 
19866831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
19870f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
198882bf6240SBarry Smith 
198982bf6240SBarry Smith   if (snes->setupcalled) {
1990e1311b90SBarry Smith     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
199182bf6240SBarry Smith     snes->data = 0;
199282bf6240SBarry Smith   }
19937d1a2b2bSBarry Smith 
19949b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
199582bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
199682bf6240SBarry Smith 
1997454a90a3SBarry Smith   ierr =  FListFind(snes->comm,SNESList,type,(int (**)(void *)) &r);CHKERRQ(ierr);
199882bf6240SBarry Smith 
1999*29bbc08cSBarry Smith   if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type);
20001548b14aSBarry Smith 
2001606d414cSSatish Balay   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
200282bf6240SBarry Smith   snes->data = 0;
20033a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
200482bf6240SBarry Smith 
2005454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
200682bf6240SBarry Smith   snes->set_method_called = 1;
200782bf6240SBarry Smith 
20083a40ed3dSBarry Smith   PetscFunctionReturn(0);
20099b94acceSBarry Smith }
20109b94acceSBarry Smith 
2011a847f771SSatish Balay 
20129b94acceSBarry Smith /* --------------------------------------------------------------------- */
20135615d1e5SSatish Balay #undef __FUNC__
2014b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESRegisterDestroy"
20159b94acceSBarry Smith /*@C
20169b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
2017f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
20189b94acceSBarry Smith 
2019fee21e36SBarry Smith    Not Collective
2020fee21e36SBarry Smith 
202136851e7fSLois Curfman McInnes    Level: advanced
202236851e7fSLois Curfman McInnes 
20239b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
20249b94acceSBarry Smith 
20259b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
20269b94acceSBarry Smith @*/
2027cf256101SBarry Smith int SNESRegisterDestroy(void)
20289b94acceSBarry Smith {
202982bf6240SBarry Smith   int ierr;
203082bf6240SBarry Smith 
20313a40ed3dSBarry Smith   PetscFunctionBegin;
203282bf6240SBarry Smith   if (SNESList) {
20331d1367b7SBarry Smith     ierr = FListDestroy(&SNESList);CHKERRQ(ierr);
203482bf6240SBarry Smith     SNESList = 0;
20359b94acceSBarry Smith   }
20364c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
20373a40ed3dSBarry Smith   PetscFunctionReturn(0);
20389b94acceSBarry Smith }
20399b94acceSBarry Smith 
20405615d1e5SSatish Balay #undef __FUNC__
2041b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetType"
20429b94acceSBarry Smith /*@C
20439a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
20449b94acceSBarry Smith 
2045c7afd0dbSLois Curfman McInnes    Not Collective
2046c7afd0dbSLois Curfman McInnes 
20479b94acceSBarry Smith    Input Parameter:
20484b0e389bSBarry Smith .  snes - nonlinear solver context
20499b94acceSBarry Smith 
20509b94acceSBarry Smith    Output Parameter:
2051454a90a3SBarry Smith .  type - SNES method (a charactor string)
20529b94acceSBarry Smith 
205336851e7fSLois Curfman McInnes    Level: intermediate
205436851e7fSLois Curfman McInnes 
2055454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
20569b94acceSBarry Smith @*/
2057454a90a3SBarry Smith int SNESGetType(SNES snes,SNESType *type)
20589b94acceSBarry Smith {
20593a40ed3dSBarry Smith   PetscFunctionBegin;
2060184914b5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2061454a90a3SBarry Smith   *type = snes->type_name;
20623a40ed3dSBarry Smith   PetscFunctionReturn(0);
20639b94acceSBarry Smith }
20649b94acceSBarry Smith 
20655615d1e5SSatish Balay #undef __FUNC__
2066b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetSolution"
20679b94acceSBarry Smith /*@C
20689b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
20699b94acceSBarry Smith    stored.
20709b94acceSBarry Smith 
2071c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2072c7afd0dbSLois Curfman McInnes 
20739b94acceSBarry Smith    Input Parameter:
20749b94acceSBarry Smith .  snes - the SNES context
20759b94acceSBarry Smith 
20769b94acceSBarry Smith    Output Parameter:
20779b94acceSBarry Smith .  x - the solution
20789b94acceSBarry Smith 
207936851e7fSLois Curfman McInnes    Level: advanced
208036851e7fSLois Curfman McInnes 
20819b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
20829b94acceSBarry Smith 
2083a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate()
20849b94acceSBarry Smith @*/
20859b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
20869b94acceSBarry Smith {
20873a40ed3dSBarry Smith   PetscFunctionBegin;
208877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
20899b94acceSBarry Smith   *x = snes->vec_sol_always;
20903a40ed3dSBarry Smith   PetscFunctionReturn(0);
20919b94acceSBarry Smith }
20929b94acceSBarry Smith 
20935615d1e5SSatish Balay #undef __FUNC__
2094b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetSolutionUpdate"
20959b94acceSBarry Smith /*@C
20969b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
20979b94acceSBarry Smith    stored.
20989b94acceSBarry Smith 
2099c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2100c7afd0dbSLois Curfman McInnes 
21019b94acceSBarry Smith    Input Parameter:
21029b94acceSBarry Smith .  snes - the SNES context
21039b94acceSBarry Smith 
21049b94acceSBarry Smith    Output Parameter:
21059b94acceSBarry Smith .  x - the solution update
21069b94acceSBarry Smith 
210736851e7fSLois Curfman McInnes    Level: advanced
210836851e7fSLois Curfman McInnes 
21099b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
21109b94acceSBarry Smith 
21119b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
21129b94acceSBarry Smith @*/
21139b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
21149b94acceSBarry Smith {
21153a40ed3dSBarry Smith   PetscFunctionBegin;
211677c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
21179b94acceSBarry Smith   *x = snes->vec_sol_update_always;
21183a40ed3dSBarry Smith   PetscFunctionReturn(0);
21199b94acceSBarry Smith }
21209b94acceSBarry Smith 
21215615d1e5SSatish Balay #undef __FUNC__
2122b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetFunction"
21239b94acceSBarry Smith /*@C
21243638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
21259b94acceSBarry Smith 
2126c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2127c7afd0dbSLois Curfman McInnes 
21289b94acceSBarry Smith    Input Parameter:
21299b94acceSBarry Smith .  snes - the SNES context
21309b94acceSBarry Smith 
21319b94acceSBarry Smith    Output Parameter:
21327bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
213300036973SBarry Smith .  ctx - the function context (or PETSC_NULL)
213400036973SBarry Smith -  func - the function (or PETSC_NULL)
21359b94acceSBarry Smith 
21369b94acceSBarry Smith    Notes:
21379b94acceSBarry Smith    SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only
21389b94acceSBarry Smith    Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are
21399b94acceSBarry Smith    SNESGetMinimizationFunction() and SNESGetGradient();
21409b94acceSBarry Smith 
214136851e7fSLois Curfman McInnes    Level: advanced
214236851e7fSLois Curfman McInnes 
2143a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
21449b94acceSBarry Smith 
214531615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(),
214631615d04SLois Curfman McInnes           SNESGetGradient()
21477bf4e008SBarry Smith 
21489b94acceSBarry Smith @*/
214900036973SBarry Smith int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*))
21509b94acceSBarry Smith {
21513a40ed3dSBarry Smith   PetscFunctionBegin;
215277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2153a8c6a408SBarry Smith   if (snes->method_class != SNES_NONLINEAR_EQUATIONS) {
2154*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_NONLINEAR_EQUATIONS only");
2155a8c6a408SBarry Smith   }
21567bf4e008SBarry Smith   if (r)    *r    = snes->vec_func_always;
21577bf4e008SBarry Smith   if (ctx)  *ctx  = snes->funP;
215800036973SBarry Smith   if (func) *func = snes->computefunction;
21593a40ed3dSBarry Smith   PetscFunctionReturn(0);
21609b94acceSBarry Smith }
21619b94acceSBarry Smith 
21625615d1e5SSatish Balay #undef __FUNC__
2163b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetGradient"
21649b94acceSBarry Smith /*@C
21653638b69dSLois Curfman McInnes    SNESGetGradient - Returns the vector where the gradient is stored.
21669b94acceSBarry Smith 
2167c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
2168c7afd0dbSLois Curfman McInnes 
21699b94acceSBarry Smith    Input Parameter:
21709b94acceSBarry Smith .  snes - the SNES context
21719b94acceSBarry Smith 
21729b94acceSBarry Smith    Output Parameter:
21737bf4e008SBarry Smith +  r - the gradient (or PETSC_NULL)
21747bf4e008SBarry Smith -  ctx - the gradient context (or PETSC_NULL)
21759b94acceSBarry Smith 
21769b94acceSBarry Smith    Notes:
21779b94acceSBarry Smith    SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods
21789b94acceSBarry Smith    only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
21799b94acceSBarry Smith    SNESGetFunction().
21809b94acceSBarry Smith 
218136851e7fSLois Curfman McInnes    Level: advanced
218236851e7fSLois Curfman McInnes 
21839b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient
21849b94acceSBarry Smith 
21857bf4e008SBarry Smith .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction(),
21867bf4e008SBarry Smith           SNESSetGradient(), SNESSetFunction()
21877bf4e008SBarry Smith 
21889b94acceSBarry Smith @*/
21897bf4e008SBarry Smith int SNESGetGradient(SNES snes,Vec *r,void **ctx)
21909b94acceSBarry Smith {
21913a40ed3dSBarry Smith   PetscFunctionBegin;
219277c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
21933a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
2194*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
21953a40ed3dSBarry Smith   }
21967bf4e008SBarry Smith   if (r)   *r = snes->vec_func_always;
21977bf4e008SBarry Smith   if (ctx) *ctx = snes->funP;
21983a40ed3dSBarry Smith   PetscFunctionReturn(0);
21999b94acceSBarry Smith }
22009b94acceSBarry Smith 
22015615d1e5SSatish Balay #undef __FUNC__
2202b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetMinimizationFunction"
22037bf4e008SBarry Smith /*@C
22049b94acceSBarry Smith    SNESGetMinimizationFunction - Returns the scalar function value for
22059b94acceSBarry Smith    unconstrained minimization problems.
22069b94acceSBarry Smith 
2207c7afd0dbSLois Curfman McInnes    Not Collective
2208c7afd0dbSLois Curfman McInnes 
22099b94acceSBarry Smith    Input Parameter:
22109b94acceSBarry Smith .  snes - the SNES context
22119b94acceSBarry Smith 
22129b94acceSBarry Smith    Output Parameter:
22137bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
22147bf4e008SBarry Smith -  ctx - the function context (or PETSC_NULL)
22159b94acceSBarry Smith 
22169b94acceSBarry Smith    Notes:
22179b94acceSBarry Smith    SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION
22189b94acceSBarry Smith    methods only.  An analogous routine for SNES_NONLINEAR_EQUATIONS methods is
22199b94acceSBarry Smith    SNESGetFunction().
22209b94acceSBarry Smith 
222136851e7fSLois Curfman McInnes    Level: advanced
222236851e7fSLois Curfman McInnes 
22239b94acceSBarry Smith .keywords: SNES, nonlinear, get, function
22249b94acceSBarry Smith 
22257bf4e008SBarry Smith .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction(), SNESSetFunction()
22267bf4e008SBarry Smith 
22279b94acceSBarry Smith @*/
2228329f5518SBarry Smith int SNESGetMinimizationFunction(SNES snes,PetscReal *r,void **ctx)
22299b94acceSBarry Smith {
22303a40ed3dSBarry Smith   PetscFunctionBegin;
223177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
223274679c65SBarry Smith   PetscValidScalarPointer(r);
22333a40ed3dSBarry Smith   if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) {
2234*29bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_WRONG,"For SNES_UNCONSTRAINED_MINIMIZATION only");
22353a40ed3dSBarry Smith   }
22367bf4e008SBarry Smith   if (r)   *r = snes->fc;
22377bf4e008SBarry Smith   if (ctx) *ctx = snes->umfunP;
22383a40ed3dSBarry Smith   PetscFunctionReturn(0);
22399b94acceSBarry Smith }
22409b94acceSBarry Smith 
22415615d1e5SSatish Balay #undef __FUNC__
2242b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESSetOptionsPrefix"
22433c7409f5SSatish Balay /*@C
22443c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2245d850072dSLois Curfman McInnes    SNES options in the database.
22463c7409f5SSatish Balay 
2247fee21e36SBarry Smith    Collective on SNES
2248fee21e36SBarry Smith 
2249c7afd0dbSLois Curfman McInnes    Input Parameter:
2250c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2251c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2252c7afd0dbSLois Curfman McInnes 
2253d850072dSLois Curfman McInnes    Notes:
2254a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2255c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2256d850072dSLois Curfman McInnes 
225736851e7fSLois Curfman McInnes    Level: advanced
225836851e7fSLois Curfman McInnes 
22593c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2260a86d99e1SLois Curfman McInnes 
2261a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
22623c7409f5SSatish Balay @*/
22633c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix)
22643c7409f5SSatish Balay {
22653c7409f5SSatish Balay   int ierr;
22663c7409f5SSatish Balay 
22673a40ed3dSBarry Smith   PetscFunctionBegin;
226877c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2269639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
22703c7409f5SSatish Balay   ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
22713a40ed3dSBarry Smith   PetscFunctionReturn(0);
22723c7409f5SSatish Balay }
22733c7409f5SSatish Balay 
22745615d1e5SSatish Balay #undef __FUNC__
2275b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESAppendOptionsPrefix"
22763c7409f5SSatish Balay /*@C
2277f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2278d850072dSLois Curfman McInnes    SNES options in the database.
22793c7409f5SSatish Balay 
2280fee21e36SBarry Smith    Collective on SNES
2281fee21e36SBarry Smith 
2282c7afd0dbSLois Curfman McInnes    Input Parameters:
2283c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2284c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2285c7afd0dbSLois Curfman McInnes 
2286d850072dSLois Curfman McInnes    Notes:
2287a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2288c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2289d850072dSLois Curfman McInnes 
229036851e7fSLois Curfman McInnes    Level: advanced
229136851e7fSLois Curfman McInnes 
22923c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2293a86d99e1SLois Curfman McInnes 
2294a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
22953c7409f5SSatish Balay @*/
22963c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix)
22973c7409f5SSatish Balay {
22983c7409f5SSatish Balay   int ierr;
22993c7409f5SSatish Balay 
23003a40ed3dSBarry Smith   PetscFunctionBegin;
230177c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2302639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
23033c7409f5SSatish Balay   ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr);
23043a40ed3dSBarry Smith   PetscFunctionReturn(0);
23053c7409f5SSatish Balay }
23063c7409f5SSatish Balay 
23075615d1e5SSatish Balay #undef __FUNC__
2308b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESGetOptionsPrefix"
23099ab63eb5SSatish Balay /*@C
23103c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
23113c7409f5SSatish Balay    SNES options in the database.
23123c7409f5SSatish Balay 
2313c7afd0dbSLois Curfman McInnes    Not Collective
2314c7afd0dbSLois Curfman McInnes 
23153c7409f5SSatish Balay    Input Parameter:
23163c7409f5SSatish Balay .  snes - the SNES context
23173c7409f5SSatish Balay 
23183c7409f5SSatish Balay    Output Parameter:
23193c7409f5SSatish Balay .  prefix - pointer to the prefix string used
23203c7409f5SSatish Balay 
23219ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
23229ab63eb5SSatish Balay    sufficient length to hold the prefix.
23239ab63eb5SSatish Balay 
232436851e7fSLois Curfman McInnes    Level: advanced
232536851e7fSLois Curfman McInnes 
23263c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2327a86d99e1SLois Curfman McInnes 
2328a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
23293c7409f5SSatish Balay @*/
23303c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix)
23313c7409f5SSatish Balay {
23323c7409f5SSatish Balay   int ierr;
23333c7409f5SSatish Balay 
23343a40ed3dSBarry Smith   PetscFunctionBegin;
233577c4ece6SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE);
2336639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
23373a40ed3dSBarry Smith   PetscFunctionReturn(0);
23383c7409f5SSatish Balay }
23393c7409f5SSatish Balay 
2340acb85ae6SSatish Balay /*MC
2341f1af5d2fSBarry Smith    SNESRegisterDynamic - Adds a method to the nonlinear solver package.
23429b94acceSBarry Smith 
2343b2002411SLois Curfman McInnes    Synopsis:
2344f1af5d2fSBarry Smith    SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES))
23459b94acceSBarry Smith 
23468d76a1e5SLois Curfman McInnes    Not collective
23478d76a1e5SLois Curfman McInnes 
2348b2002411SLois Curfman McInnes    Input Parameters:
2349c7afd0dbSLois Curfman McInnes +  name_solver - name of a new user-defined solver
2350b2002411SLois Curfman McInnes .  path - path (either absolute or relative) the library containing this solver
2351b2002411SLois Curfman McInnes .  name_create - name of routine to create method context
2352c7afd0dbSLois Curfman McInnes -  routine_create - routine to create method context
23539b94acceSBarry Smith 
2354b2002411SLois Curfman McInnes    Notes:
2355f1af5d2fSBarry Smith    SNESRegisterDynamic() may be called multiple times to add several user-defined solvers.
23563a40ed3dSBarry Smith 
2357b2002411SLois Curfman McInnes    If dynamic libraries are used, then the fourth input argument (routine_create)
2358b2002411SLois Curfman McInnes    is ignored.
2359b2002411SLois Curfman McInnes 
23605ba88a07SLois Curfman McInnes    Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LDIR}, ${BOPT},
23615ba88a07SLois Curfman McInnes    and others of the form ${any_environmental_variable} occuring in pathname will be
23625ba88a07SLois Curfman McInnes    replaced with appropriate values.
23635ba88a07SLois Curfman McInnes 
2364b2002411SLois Curfman McInnes    Sample usage:
236569225d78SLois Curfman McInnes .vb
2366f1af5d2fSBarry Smith    SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a,
2367b2002411SLois Curfman McInnes                 "MySolverCreate",MySolverCreate);
236869225d78SLois Curfman McInnes .ve
2369b2002411SLois Curfman McInnes 
2370b2002411SLois Curfman McInnes    Then, your solver can be chosen with the procedural interface via
2371b2002411SLois Curfman McInnes $     SNESSetType(snes,"my_solver")
2372b2002411SLois Curfman McInnes    or at runtime via the option
2373b2002411SLois Curfman McInnes $     -snes_type my_solver
2374b2002411SLois Curfman McInnes 
237536851e7fSLois Curfman McInnes    Level: advanced
237636851e7fSLois Curfman McInnes 
2377b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register
2378b2002411SLois Curfman McInnes 
2379b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy()
2380b2002411SLois Curfman McInnes M*/
2381b2002411SLois Curfman McInnes 
2382b2002411SLois Curfman McInnes #undef __FUNC__
2383b2863d3aSBarry Smith #define __FUNC__ /*<a name=""></a>*/"SNESRegister"
2384f1af5d2fSBarry Smith int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES))
2385b2002411SLois Curfman McInnes {
2386b2002411SLois Curfman McInnes   char fullname[256];
2387b2002411SLois Curfman McInnes   int  ierr;
2388b2002411SLois Curfman McInnes 
2389b2002411SLois Curfman McInnes   PetscFunctionBegin;
2390f1af5d2fSBarry Smith   ierr = FListConcat(path,name,fullname);CHKERRQ(ierr);
2391f1af5d2fSBarry Smith   ierr = FListAdd(&SNESList,sname,fullname,(int (*)(void*))function);CHKERRQ(ierr);
2392b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2393b2002411SLois Curfman McInnes }
2394