xref: /petsc/src/snes/interface/snes.c (revision 32077d6df18512e94702c86e7d562481ed0cebd0)
19b94acceSBarry Smith 
2e090d566SSatish Balay #include "src/snes/snesimpl.h"      /*I "petscsnes.h"  I*/
39b94acceSBarry Smith 
44c49b128SBarry Smith PetscTruth SNESRegisterAllCalled = PETSC_FALSE;
58ba1e511SMatthew Knepley PetscFList SNESList              = PETSC_NULL;
68ba1e511SMatthew Knepley 
78ba1e511SMatthew Knepley /* Logging support */
843c77886SBarry Smith int SNES_COOKIE = 0;
943c77886SBarry Smith int SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0;
10a09944afSBarry Smith 
11a09944afSBarry Smith #undef __FUNCT__
124a2ae208SSatish Balay #define __FUNCT__ "SNESView"
137e2c5f70SBarry Smith /*@C
149b94acceSBarry Smith    SNESView - Prints the SNES data structure.
159b94acceSBarry Smith 
164c49b128SBarry Smith    Collective on SNES
17fee21e36SBarry Smith 
18c7afd0dbSLois Curfman McInnes    Input Parameters:
19c7afd0dbSLois Curfman McInnes +  SNES - the SNES context
20c7afd0dbSLois Curfman McInnes -  viewer - visualization context
21c7afd0dbSLois Curfman McInnes 
229b94acceSBarry Smith    Options Database Key:
23c8a8ba5cSLois Curfman McInnes .  -snes_view - Calls SNESView() at end of SNESSolve()
249b94acceSBarry Smith 
259b94acceSBarry Smith    Notes:
269b94acceSBarry Smith    The available visualization contexts include
27b0a32e0cSBarry Smith +     PETSC_VIEWER_STDOUT_SELF - standard output (default)
28b0a32e0cSBarry Smith -     PETSC_VIEWER_STDOUT_WORLD - synchronized standard
29c8a8ba5cSLois Curfman McInnes          output where only the first processor opens
30c8a8ba5cSLois Curfman McInnes          the file.  All other processors send their
31c8a8ba5cSLois Curfman McInnes          data to the first processor to print.
329b94acceSBarry Smith 
333e081fefSLois Curfman McInnes    The user can open an alternative visualization context with
34b0a32e0cSBarry Smith    PetscViewerASCIIOpen() - output to a specified file.
359b94acceSBarry Smith 
3636851e7fSLois Curfman McInnes    Level: beginner
3736851e7fSLois Curfman McInnes 
389b94acceSBarry Smith .keywords: SNES, view
399b94acceSBarry Smith 
40b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen()
419b94acceSBarry Smith @*/
42b0a32e0cSBarry Smith int SNESView(SNES snes,PetscViewer viewer)
439b94acceSBarry Smith {
449b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
459b94acceSBarry Smith   int                 ierr;
4694b7f48cSBarry Smith   KSP                 ksp;
47454a90a3SBarry Smith   char                *type;
48*32077d6dSBarry Smith   PetscTruth          iascii,isstring;
499b94acceSBarry Smith 
503a40ed3dSBarry Smith   PetscFunctionBegin;
514482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
52b0a32e0cSBarry Smith   if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm);
534482741eSBarry Smith   PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2);
54c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,viewer,2);
5574679c65SBarry Smith 
56*32077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
57b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
58*32077d6dSBarry Smith   if (iascii) {
593a7fca6bSBarry Smith     if (snes->prefix) {
603a7fca6bSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr);
613a7fca6bSBarry Smith     } else {
62b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr);
633a7fca6bSBarry Smith     }
64454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
65454a90a3SBarry Smith     if (type) {
66b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: %s\n",type);CHKERRQ(ierr);
67184914b5SBarry Smith     } else {
68b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer,"  type: not set yet\n");CHKERRQ(ierr);
69184914b5SBarry Smith     }
700ef38995SBarry Smith     if (snes->view) {
71b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
720ef38995SBarry Smith       ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr);
73b0a32e0cSBarry Smith       ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
740ef38995SBarry Smith     }
75b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr);
7677d8c4bbSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  tolerances: relative=%g, absolute=%g, solution=%g\n",
7777d8c4bbSBarry Smith                  snes->rtol,snes->atol,snes->xtol);CHKERRQ(ierr);
78b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr);
79b0a32e0cSBarry Smith     ierr = PetscViewerASCIIPrintf(viewer,"  total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr);
809b94acceSBarry Smith     if (snes->ksp_ewconv) {
819b94acceSBarry Smith       kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
829b94acceSBarry Smith       if (kctx) {
83b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"  Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr);
84b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr);
85b0a32e0cSBarry Smith         ierr = PetscViewerASCIIPrintf(viewer,"    gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr);
869b94acceSBarry Smith       }
879b94acceSBarry Smith     }
880f5bd95cSBarry Smith   } else if (isstring) {
89454a90a3SBarry Smith     ierr = SNESGetType(snes,&type);CHKERRQ(ierr);
90b0a32e0cSBarry Smith     ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr);
9119bcc07fSBarry Smith   }
9294b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
93b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
9494b7f48cSBarry Smith   ierr = KSPView(ksp,viewer);CHKERRQ(ierr);
95b0a32e0cSBarry Smith   ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
963a40ed3dSBarry Smith   PetscFunctionReturn(0);
979b94acceSBarry Smith }
989b94acceSBarry Smith 
9976b2cf59SMatthew Knepley /*
10076b2cf59SMatthew Knepley   We retain a list of functions that also take SNES command
10176b2cf59SMatthew Knepley   line options. These are called at the end SNESSetFromOptions()
10276b2cf59SMatthew Knepley */
10376b2cf59SMatthew Knepley #define MAXSETFROMOPTIONS 5
10476b2cf59SMatthew Knepley static int numberofsetfromoptions;
10576b2cf59SMatthew Knepley static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES);
10676b2cf59SMatthew Knepley 
107e74ef692SMatthew Knepley #undef __FUNCT__
108e74ef692SMatthew Knepley #define __FUNCT__ "SNESAddOptionsChecker"
109ac226902SBarry Smith /*@C
11076b2cf59SMatthew Knepley   SNESAddOptionsChecker - Adds an additional function to check for SNES options.
11176b2cf59SMatthew Knepley 
11276b2cf59SMatthew Knepley   Not Collective
11376b2cf59SMatthew Knepley 
11476b2cf59SMatthew Knepley   Input Parameter:
11576b2cf59SMatthew Knepley . snescheck - function that checks for options
11676b2cf59SMatthew Knepley 
11776b2cf59SMatthew Knepley   Level: developer
11876b2cf59SMatthew Knepley 
11976b2cf59SMatthew Knepley .seealso: SNESSetFromOptions()
12076b2cf59SMatthew Knepley @*/
12176b2cf59SMatthew Knepley int SNESAddOptionsChecker(int (*snescheck)(SNES))
12276b2cf59SMatthew Knepley {
12376b2cf59SMatthew Knepley   PetscFunctionBegin;
12476b2cf59SMatthew Knepley   if (numberofsetfromoptions >= MAXSETFROMOPTIONS) {
12576b2cf59SMatthew Knepley     SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %d allowed", MAXSETFROMOPTIONS);
12676b2cf59SMatthew Knepley   }
12776b2cf59SMatthew Knepley 
12876b2cf59SMatthew Knepley   othersetfromoptions[numberofsetfromoptions++] = snescheck;
12976b2cf59SMatthew Knepley   PetscFunctionReturn(0);
13076b2cf59SMatthew Knepley }
13176b2cf59SMatthew Knepley 
1324a2ae208SSatish Balay #undef __FUNCT__
1334a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions"
1349b94acceSBarry Smith /*@
13594b7f48cSBarry Smith    SNESSetFromOptions - Sets various SNES and KSP parameters from user options.
1369b94acceSBarry Smith 
137c7afd0dbSLois Curfman McInnes    Collective on SNES
138c7afd0dbSLois Curfman McInnes 
1399b94acceSBarry Smith    Input Parameter:
1409b94acceSBarry Smith .  snes - the SNES context
1419b94acceSBarry Smith 
14236851e7fSLois Curfman McInnes    Options Database Keys:
1436831982aSBarry Smith +  -snes_type <type> - ls, tr, umls, umtr, test
14482738288SBarry Smith .  -snes_stol - convergence tolerance in terms of the norm
14582738288SBarry Smith                 of the change in the solution between steps
146b39c3a46SLois Curfman McInnes .  -snes_atol <atol> - absolute tolerance of residual norm
147b39c3a46SLois Curfman McInnes .  -snes_rtol <rtol> - relative decrease in tolerance norm from initial
148b39c3a46SLois Curfman McInnes .  -snes_max_it <max_it> - maximum number of iterations
149b39c3a46SLois Curfman McInnes .  -snes_max_funcs <max_funcs> - maximum number of function evaluations
15050ffb88aSMatthew Knepley .  -snes_max_fail <max_fail> - maximum number of failures
151b39c3a46SLois Curfman McInnes .  -snes_trtol <trtol> - trust region tolerance
15282738288SBarry Smith .  -snes_no_convergence_test - skip convergence test in nonlinear or minimization
15382738288SBarry Smith                                solver; hence iterations will continue until max_it
1541fbbfb26SLois Curfman McInnes                                or some other criterion is reached. Saves expense
15582738288SBarry Smith                                of convergence test
15682738288SBarry Smith .  -snes_monitor - prints residual norm at each iteration
157d132466eSBarry Smith .  -snes_vecmonitor - plots solution at each iteration
158d132466eSBarry Smith .  -snes_vecmonitor_update - plots update to solution at each iteration
15982738288SBarry Smith .  -snes_xmonitor - plots residual norm at each iteration
160e24b481bSBarry Smith .  -snes_fd - use finite differences to compute Jacobian; very slow, only for testing
1615968eb51SBarry Smith .  -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration
1625968eb51SBarry Smith -  -snes_print_converged_reason - print the reason for convergence/divergence after each solve
16382738288SBarry Smith 
16482738288SBarry Smith     Options Database for Eisenstat-Walker method:
1654b27c08aSLois Curfman McInnes +  -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence
1664b27c08aSLois Curfman McInnes .  -snes_ksp_ew_version ver - version of  Eisenstat-Walker method
16736851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0
16836851e7fSLois Curfman McInnes .  -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax
16936851e7fSLois Curfman McInnes .  -snes_ksp_ew_gamma <gamma> - Sets gamma
17036851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha <alpha> - Sets alpha
17136851e7fSLois Curfman McInnes .  -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2
17236851e7fSLois Curfman McInnes -  -snes_ksp_ew_threshold <threshold> - Sets threshold
17382738288SBarry Smith 
17411ca99fdSLois Curfman McInnes    Notes:
17511ca99fdSLois Curfman McInnes    To see all options, run your program with the -help option or consult
17611ca99fdSLois Curfman McInnes    the users manual.
17783e2fdc7SBarry Smith 
17836851e7fSLois Curfman McInnes    Level: beginner
17936851e7fSLois Curfman McInnes 
1809b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database
1819b94acceSBarry Smith 
18269ed319cSSatish Balay .seealso: SNESSetOptionsPrefix()
1839b94acceSBarry Smith @*/
1849b94acceSBarry Smith int SNESSetFromOptions(SNES snes)
1859b94acceSBarry Smith {
18694b7f48cSBarry Smith   KSP                 ksp;
187186905e3SBarry Smith   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
188f1af5d2fSBarry Smith   PetscTruth          flg;
18976b2cf59SMatthew Knepley   int                 ierr, i;
1902fc52814SBarry Smith   const char          *deft;
1912fc52814SBarry Smith   char                type[256];
1929b94acceSBarry Smith 
1933a40ed3dSBarry Smith   PetscFunctionBegin;
1944482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
195ca161407SBarry Smith 
196b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
197186905e3SBarry Smith     if (snes->type_name) {
198186905e3SBarry Smith       deft = snes->type_name;
199186905e3SBarry Smith     } else {
2004b27c08aSLois Curfman McInnes       deft = SNESLS;
201d64ed03dSBarry Smith     }
2024bbc92c1SBarry Smith 
203186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
204b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
205d64ed03dSBarry Smith     if (flg) {
206186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
207186905e3SBarry Smith     } else if (!snes->type_name) {
208186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
209d64ed03dSBarry Smith     }
210909c8a9fSBarry Smith     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
21193c39befSBarry Smith 
21287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
21387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr);
214186905e3SBarry Smith 
21587828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
216b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
217b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
21850ffb88aSMatthew Knepley     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
2195968eb51SBarry Smith     ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr);
2205968eb51SBarry Smith     if (flg) {
2215968eb51SBarry Smith       snes->printreason = PETSC_TRUE;
2225968eb51SBarry Smith     }
223186905e3SBarry Smith 
224b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr);
225186905e3SBarry Smith 
226b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
22787828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
22887828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
22987828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
23087828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
23187828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
23287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
233186905e3SBarry Smith 
234b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
23593c39befSBarry Smith     if (flg) {snes->converged = 0;}
236b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
2375cd90555SBarry Smith     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
238b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
239b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
2403a7fca6bSBarry Smith     ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr);
2413a7fca6bSBarry Smith     if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);}
242b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
243b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
244b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
245b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
246b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
2477c922b88SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
2485ed2d596SBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr);
2495ed2d596SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);}
250b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
251186905e3SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
252e24b481bSBarry Smith 
253b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
2544b27c08aSLois Curfman McInnes     if (flg) {
255186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
256b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
2579b94acceSBarry Smith     }
258639f9d9dSBarry Smith 
25976b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
26076b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
26176b2cf59SMatthew Knepley     }
26276b2cf59SMatthew Knepley 
263186905e3SBarry Smith     if (snes->setfromoptions) {
264186905e3SBarry Smith       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
265639f9d9dSBarry Smith     }
266639f9d9dSBarry Smith 
267b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2684bbc92c1SBarry Smith 
26994b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
27094b7f48cSBarry Smith   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
27193993e2dSLois Curfman McInnes 
2723a40ed3dSBarry Smith   PetscFunctionReturn(0);
2739b94acceSBarry Smith }
2749b94acceSBarry Smith 
275a847f771SSatish Balay 
2764a2ae208SSatish Balay #undef __FUNCT__
2774a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
2789b94acceSBarry Smith /*@
2799b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2809b94acceSBarry Smith    the nonlinear solvers.
2819b94acceSBarry Smith 
282fee21e36SBarry Smith    Collective on SNES
283fee21e36SBarry Smith 
284c7afd0dbSLois Curfman McInnes    Input Parameters:
285c7afd0dbSLois Curfman McInnes +  snes - the SNES context
286c7afd0dbSLois Curfman McInnes -  usrP - optional user context
287c7afd0dbSLois Curfman McInnes 
28836851e7fSLois Curfman McInnes    Level: intermediate
28936851e7fSLois Curfman McInnes 
2909b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2919b94acceSBarry Smith 
2929b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2939b94acceSBarry Smith @*/
2949b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP)
2959b94acceSBarry Smith {
2963a40ed3dSBarry Smith   PetscFunctionBegin;
2974482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2989b94acceSBarry Smith   snes->user		= usrP;
2993a40ed3dSBarry Smith   PetscFunctionReturn(0);
3009b94acceSBarry Smith }
30174679c65SBarry Smith 
3024a2ae208SSatish Balay #undef __FUNCT__
3034a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
3049b94acceSBarry Smith /*@C
3059b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3069b94acceSBarry Smith    nonlinear solvers.
3079b94acceSBarry Smith 
308c7afd0dbSLois Curfman McInnes    Not Collective
309c7afd0dbSLois Curfman McInnes 
3109b94acceSBarry Smith    Input Parameter:
3119b94acceSBarry Smith .  snes - SNES context
3129b94acceSBarry Smith 
3139b94acceSBarry Smith    Output Parameter:
3149b94acceSBarry Smith .  usrP - user context
3159b94acceSBarry Smith 
31636851e7fSLois Curfman McInnes    Level: intermediate
31736851e7fSLois Curfman McInnes 
3189b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3199b94acceSBarry Smith 
3209b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3219b94acceSBarry Smith @*/
3229b94acceSBarry Smith int SNESGetApplicationContext(SNES snes,void **usrP)
3239b94acceSBarry Smith {
3243a40ed3dSBarry Smith   PetscFunctionBegin;
3254482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3269b94acceSBarry Smith   *usrP = snes->user;
3273a40ed3dSBarry Smith   PetscFunctionReturn(0);
3289b94acceSBarry Smith }
32974679c65SBarry Smith 
3304a2ae208SSatish Balay #undef __FUNCT__
3314a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
3329b94acceSBarry Smith /*@
333c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
334c8228a4eSBarry Smith    at this time.
3359b94acceSBarry Smith 
336c7afd0dbSLois Curfman McInnes    Not Collective
337c7afd0dbSLois Curfman McInnes 
3389b94acceSBarry Smith    Input Parameter:
3399b94acceSBarry Smith .  snes - SNES context
3409b94acceSBarry Smith 
3419b94acceSBarry Smith    Output Parameter:
3429b94acceSBarry Smith .  iter - iteration number
3439b94acceSBarry Smith 
344c8228a4eSBarry Smith    Notes:
345c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
346c8228a4eSBarry Smith 
347c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
34808405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
34908405cd6SLois Curfman McInnes .vb
35008405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
35108405cd6SLois Curfman McInnes       if (!(it % 2)) {
35208405cd6SLois Curfman McInnes         [compute Jacobian here]
35308405cd6SLois Curfman McInnes       }
35408405cd6SLois Curfman McInnes .ve
355c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
35608405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
357c8228a4eSBarry Smith 
35836851e7fSLois Curfman McInnes    Level: intermediate
35936851e7fSLois Curfman McInnes 
3609b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3619b94acceSBarry Smith @*/
3629b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter)
3639b94acceSBarry Smith {
3643a40ed3dSBarry Smith   PetscFunctionBegin;
3654482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3664482741eSBarry Smith   PetscValidIntPointer(iter,2);
3679b94acceSBarry Smith   *iter = snes->iter;
3683a40ed3dSBarry Smith   PetscFunctionReturn(0);
3699b94acceSBarry Smith }
37074679c65SBarry Smith 
3714a2ae208SSatish Balay #undef __FUNCT__
3724a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
3739b94acceSBarry Smith /*@
3749b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3759b94acceSBarry Smith    with SNESSSetFunction().
3769b94acceSBarry Smith 
377c7afd0dbSLois Curfman McInnes    Collective on SNES
378c7afd0dbSLois Curfman McInnes 
3799b94acceSBarry Smith    Input Parameter:
3809b94acceSBarry Smith .  snes - SNES context
3819b94acceSBarry Smith 
3829b94acceSBarry Smith    Output Parameter:
3839b94acceSBarry Smith .  fnorm - 2-norm of function
3849b94acceSBarry Smith 
38536851e7fSLois Curfman McInnes    Level: intermediate
38636851e7fSLois Curfman McInnes 
3879b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
388a86d99e1SLois Curfman McInnes 
38908405cd6SLois Curfman McInnes .seealso: SNESGetFunction()
3909b94acceSBarry Smith @*/
39187828ca2SBarry Smith int SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
3929b94acceSBarry Smith {
3933a40ed3dSBarry Smith   PetscFunctionBegin;
3944482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3954482741eSBarry Smith   PetscValidScalarPointer(fnorm,2);
3969b94acceSBarry Smith   *fnorm = snes->norm;
3973a40ed3dSBarry Smith   PetscFunctionReturn(0);
3989b94acceSBarry Smith }
39974679c65SBarry Smith 
4004a2ae208SSatish Balay #undef __FUNCT__
4014a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
4029b94acceSBarry Smith /*@
4039b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4049b94acceSBarry Smith    attempted by the nonlinear solver.
4059b94acceSBarry Smith 
406c7afd0dbSLois Curfman McInnes    Not Collective
407c7afd0dbSLois Curfman McInnes 
4089b94acceSBarry Smith    Input Parameter:
4099b94acceSBarry Smith .  snes - SNES context
4109b94acceSBarry Smith 
4119b94acceSBarry Smith    Output Parameter:
4129b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4139b94acceSBarry Smith 
414c96a6f78SLois Curfman McInnes    Notes:
415c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
416c96a6f78SLois Curfman McInnes 
41736851e7fSLois Curfman McInnes    Level: intermediate
41836851e7fSLois Curfman McInnes 
4199b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4209b94acceSBarry Smith @*/
4219b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4229b94acceSBarry Smith {
4233a40ed3dSBarry Smith   PetscFunctionBegin;
4244482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4254482741eSBarry Smith   PetscValidIntPointer(nfails,2);
42650ffb88aSMatthew Knepley   *nfails = snes->numFailures;
42750ffb88aSMatthew Knepley   PetscFunctionReturn(0);
42850ffb88aSMatthew Knepley }
42950ffb88aSMatthew Knepley 
43050ffb88aSMatthew Knepley #undef __FUNCT__
43150ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps"
43250ffb88aSMatthew Knepley /*@
43350ffb88aSMatthew Knepley    SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps
43450ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
43550ffb88aSMatthew Knepley 
43650ffb88aSMatthew Knepley    Not Collective
43750ffb88aSMatthew Knepley 
43850ffb88aSMatthew Knepley    Input Parameters:
43950ffb88aSMatthew Knepley +  snes     - SNES context
44050ffb88aSMatthew Knepley -  maxFails - maximum of unsuccessful steps
44150ffb88aSMatthew Knepley 
44250ffb88aSMatthew Knepley    Level: intermediate
44350ffb88aSMatthew Knepley 
44450ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
44550ffb88aSMatthew Knepley @*/
44650ffb88aSMatthew Knepley int SNESSetMaximumUnsuccessfulSteps(SNES snes, int maxFails)
44750ffb88aSMatthew Knepley {
44850ffb88aSMatthew Knepley   PetscFunctionBegin;
4494482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
45050ffb88aSMatthew Knepley   snes->maxFailures = maxFails;
45150ffb88aSMatthew Knepley   PetscFunctionReturn(0);
45250ffb88aSMatthew Knepley }
45350ffb88aSMatthew Knepley 
45450ffb88aSMatthew Knepley #undef __FUNCT__
45550ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps"
45650ffb88aSMatthew Knepley /*@
45750ffb88aSMatthew Knepley    SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps
45850ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
45950ffb88aSMatthew Knepley 
46050ffb88aSMatthew Knepley    Not Collective
46150ffb88aSMatthew Knepley 
46250ffb88aSMatthew Knepley    Input Parameter:
46350ffb88aSMatthew Knepley .  snes     - SNES context
46450ffb88aSMatthew Knepley 
46550ffb88aSMatthew Knepley    Output Parameter:
46650ffb88aSMatthew Knepley .  maxFails - maximum of unsuccessful steps
46750ffb88aSMatthew Knepley 
46850ffb88aSMatthew Knepley    Level: intermediate
46950ffb88aSMatthew Knepley 
47050ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
47150ffb88aSMatthew Knepley @*/
47250ffb88aSMatthew Knepley int SNESGetMaximumUnsuccessfulSteps(SNES snes, int *maxFails)
47350ffb88aSMatthew Knepley {
47450ffb88aSMatthew Knepley   PetscFunctionBegin;
4754482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4764482741eSBarry Smith   PetscValidIntPointer(maxFails,2);
47750ffb88aSMatthew Knepley   *maxFails = snes->maxFailures;
4783a40ed3dSBarry Smith   PetscFunctionReturn(0);
4799b94acceSBarry Smith }
480a847f771SSatish Balay 
4814a2ae208SSatish Balay #undef __FUNCT__
4824a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations"
483c96a6f78SLois Curfman McInnes /*@
484c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
485c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
486c96a6f78SLois Curfman McInnes 
487c7afd0dbSLois Curfman McInnes    Not Collective
488c7afd0dbSLois Curfman McInnes 
489c96a6f78SLois Curfman McInnes    Input Parameter:
490c96a6f78SLois Curfman McInnes .  snes - SNES context
491c96a6f78SLois Curfman McInnes 
492c96a6f78SLois Curfman McInnes    Output Parameter:
493c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
494c96a6f78SLois Curfman McInnes 
495c96a6f78SLois Curfman McInnes    Notes:
496c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
497c96a6f78SLois Curfman McInnes 
49836851e7fSLois Curfman McInnes    Level: intermediate
49936851e7fSLois Curfman McInnes 
500c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
501c96a6f78SLois Curfman McInnes @*/
50286bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits)
503c96a6f78SLois Curfman McInnes {
5043a40ed3dSBarry Smith   PetscFunctionBegin;
5054482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5064482741eSBarry Smith   PetscValidIntPointer(lits,2);
507c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
5083a40ed3dSBarry Smith   PetscFunctionReturn(0);
509c96a6f78SLois Curfman McInnes }
510c96a6f78SLois Curfman McInnes 
5114a2ae208SSatish Balay #undef __FUNCT__
51294b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP"
5139b94acceSBarry Smith /*@C
51494b7f48cSBarry Smith    SNESGetKSP - Returns the KSP context for a SNES solver.
5159b94acceSBarry Smith 
51694b7f48cSBarry Smith    Not Collective, but if SNES object is parallel, then KSP object is parallel
517c7afd0dbSLois Curfman McInnes 
5189b94acceSBarry Smith    Input Parameter:
5199b94acceSBarry Smith .  snes - the SNES context
5209b94acceSBarry Smith 
5219b94acceSBarry Smith    Output Parameter:
52294b7f48cSBarry Smith .  ksp - the KSP context
5239b94acceSBarry Smith 
5249b94acceSBarry Smith    Notes:
52594b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
5269b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5279b94acceSBarry Smith    KSP and PC contexts as well.
5289b94acceSBarry Smith 
52936851e7fSLois Curfman McInnes    Level: beginner
53036851e7fSLois Curfman McInnes 
53194b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
5329b94acceSBarry Smith 
53394b7f48cSBarry Smith .seealso: KSPGetPC()
5349b94acceSBarry Smith @*/
53594b7f48cSBarry Smith int SNESGetKSP(SNES snes,KSP *ksp)
5369b94acceSBarry Smith {
5373a40ed3dSBarry Smith   PetscFunctionBegin;
5384482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5394482741eSBarry Smith   PetscValidPointer(ksp,2);
54094b7f48cSBarry Smith   *ksp = snes->ksp;
5413a40ed3dSBarry Smith   PetscFunctionReturn(0);
5429b94acceSBarry Smith }
54382bf6240SBarry Smith 
5444a2ae208SSatish Balay #undef __FUNCT__
5454a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
546454a90a3SBarry Smith static int SNESPublish_Petsc(PetscObject obj)
547e24b481bSBarry Smith {
548aa482453SBarry Smith #if defined(PETSC_HAVE_AMS)
549454a90a3SBarry Smith   SNES          v = (SNES) obj;
550e24b481bSBarry Smith   int          ierr;
55143d6d2cbSBarry Smith #endif
552e24b481bSBarry Smith 
553e24b481bSBarry Smith   PetscFunctionBegin;
554e24b481bSBarry Smith 
55543d6d2cbSBarry Smith #if defined(PETSC_HAVE_AMS)
556e24b481bSBarry Smith   /* if it is already published then return */
557e24b481bSBarry Smith   if (v->amem >=0) PetscFunctionReturn(0);
558e24b481bSBarry Smith 
559454a90a3SBarry Smith   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
560e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
561e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
562e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
563e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
564454a90a3SBarry Smith   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
565433994e6SBarry Smith #endif
566e24b481bSBarry Smith   PetscFunctionReturn(0);
567e24b481bSBarry Smith }
568e24b481bSBarry Smith 
5699b94acceSBarry Smith /* -----------------------------------------------------------*/
5704a2ae208SSatish Balay #undef __FUNCT__
5714a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
5729b94acceSBarry Smith /*@C
5739b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5749b94acceSBarry Smith 
575c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
576c7afd0dbSLois Curfman McInnes 
577c7afd0dbSLois Curfman McInnes    Input Parameters:
578c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
5799b94acceSBarry Smith 
5809b94acceSBarry Smith    Output Parameter:
5819b94acceSBarry Smith .  outsnes - the new SNES context
5829b94acceSBarry Smith 
583c7afd0dbSLois Curfman McInnes    Options Database Keys:
584c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
585c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
586c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
587c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
588c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
589c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
590c1f60f51SBarry Smith 
59136851e7fSLois Curfman McInnes    Level: beginner
59236851e7fSLois Curfman McInnes 
5939b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5949b94acceSBarry Smith 
5954b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES
5969b94acceSBarry Smith @*/
5974b27c08aSLois Curfman McInnes int SNESCreate(MPI_Comm comm,SNES *outsnes)
5989b94acceSBarry Smith {
5999b94acceSBarry Smith   int                 ierr;
6009b94acceSBarry Smith   SNES                snes;
6019b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
60237fcc0dbSBarry Smith 
6033a40ed3dSBarry Smith   PetscFunctionBegin;
6044482741eSBarry Smith   PetscValidPointer(outsnes,1);
6058ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
6068ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
6078ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
6088ba1e511SMatthew Knepley #endif
6098ba1e511SMatthew Knepley 
6103f1db9ecSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);
611b0a32e0cSBarry Smith   PetscLogObjectCreate(snes);
612e24b481bSBarry Smith   snes->bops->publish     = SNESPublish_Petsc;
6139b94acceSBarry Smith   snes->max_its           = 50;
6149750a799SBarry Smith   snes->max_funcs	  = 10000;
6159b94acceSBarry Smith   snes->norm		  = 0.0;
616b4874afaSBarry Smith   snes->rtol		  = 1.e-8;
617b4874afaSBarry Smith   snes->ttol              = 0.0;
618b4874afaSBarry Smith   snes->atol		  = 1.e-50;
6199b94acceSBarry Smith   snes->xtol		  = 1.e-8;
6204b27c08aSLois Curfman McInnes   snes->deltatol	  = 1.e-12;
6219b94acceSBarry Smith   snes->nfuncs            = 0;
62250ffb88aSMatthew Knepley   snes->numFailures       = 0;
62350ffb88aSMatthew Knepley   snes->maxFailures       = 1;
6247a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
625639f9d9dSBarry Smith   snes->numbermonitors    = 0;
6269b94acceSBarry Smith   snes->data              = 0;
6279b94acceSBarry Smith   snes->view              = 0;
62882bf6240SBarry Smith   snes->setupcalled       = 0;
629186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
6306f24a144SLois Curfman McInnes   snes->vwork             = 0;
6316f24a144SLois Curfman McInnes   snes->nwork             = 0;
632758f92a0SBarry Smith   snes->conv_hist_len     = 0;
633758f92a0SBarry Smith   snes->conv_hist_max     = 0;
634758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
635758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
636758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
637184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
6389b94acceSBarry Smith 
6399b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
640b0a32e0cSBarry Smith   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
641b0a32e0cSBarry Smith   PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
6429b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6439b94acceSBarry Smith   kctx->version     = 2;
6449b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6459b94acceSBarry Smith                              this was too large for some test cases */
6469b94acceSBarry Smith   kctx->rtol_last   = 0;
6479b94acceSBarry Smith   kctx->rtol_max    = .9;
6489b94acceSBarry Smith   kctx->gamma       = 1.0;
6499b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6509b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6519b94acceSBarry Smith   kctx->threshold   = .1;
6529b94acceSBarry Smith   kctx->lresid_last = 0;
6539b94acceSBarry Smith   kctx->norm_last   = 0;
6549b94acceSBarry Smith 
65594b7f48cSBarry Smith   ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr);
65694b7f48cSBarry Smith   PetscLogObjectParent(snes,snes->ksp)
6575334005bSBarry Smith 
6589b94acceSBarry Smith   *outsnes = snes;
65900036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
6603a40ed3dSBarry Smith   PetscFunctionReturn(0);
6619b94acceSBarry Smith }
6629b94acceSBarry Smith 
6634a2ae208SSatish Balay #undef __FUNCT__
6644a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
6659b94acceSBarry Smith /*@C
6669b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6679b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6689b94acceSBarry Smith    equations.
6699b94acceSBarry Smith 
670fee21e36SBarry Smith    Collective on SNES
671fee21e36SBarry Smith 
672c7afd0dbSLois Curfman McInnes    Input Parameters:
673c7afd0dbSLois Curfman McInnes +  snes - the SNES context
674c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
675c7afd0dbSLois Curfman McInnes .  r - vector to store function value
676c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
677c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6789b94acceSBarry Smith 
679c7afd0dbSLois Curfman McInnes    Calling sequence of func:
6808d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
681c7afd0dbSLois Curfman McInnes 
682313e4042SLois Curfman McInnes .  f - function vector
683c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
6849b94acceSBarry Smith 
6859b94acceSBarry Smith    Notes:
6869b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6879b94acceSBarry Smith $      f'(x) x = -f(x),
688c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
6899b94acceSBarry Smith 
69036851e7fSLois Curfman McInnes    Level: beginner
69136851e7fSLois Curfman McInnes 
6929b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6939b94acceSBarry Smith 
694a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6959b94acceSBarry Smith @*/
6965334005bSBarry Smith int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx)
6979b94acceSBarry Smith {
6983a40ed3dSBarry Smith   PetscFunctionBegin;
6994482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7004482741eSBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
701c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,r,2);
702184914b5SBarry Smith 
7039b94acceSBarry Smith   snes->computefunction     = func;
7049b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7059b94acceSBarry Smith   snes->funP                = ctx;
7063a40ed3dSBarry Smith   PetscFunctionReturn(0);
7079b94acceSBarry Smith }
7089b94acceSBarry Smith 
7093ab0aad5SBarry Smith /* --------------------------------------------------------------- */
7103ab0aad5SBarry Smith #undef __FUNCT__
7113ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs"
7123ab0aad5SBarry Smith /*@C
7133ab0aad5SBarry Smith    SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set
7143ab0aad5SBarry Smith    it assumes a zero right hand side.
7153ab0aad5SBarry Smith 
7163ab0aad5SBarry Smith    Collective on SNES
7173ab0aad5SBarry Smith 
7183ab0aad5SBarry Smith    Input Parameters:
7193ab0aad5SBarry Smith +  snes - the SNES context
7203ab0aad5SBarry Smith -  rhs - the right hand side vector or PETSC_NULL for a zero right hand side
7213ab0aad5SBarry Smith 
7223ab0aad5SBarry Smith    Level: intermediate
7233ab0aad5SBarry Smith 
7243ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side
7253ab0aad5SBarry Smith 
7263ab0aad5SBarry Smith .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
7273ab0aad5SBarry Smith @*/
7283ab0aad5SBarry Smith int SNESSetRhs(SNES snes,Vec rhs)
7293ab0aad5SBarry Smith {
7303ab0aad5SBarry Smith   int ierr;
7313ab0aad5SBarry Smith 
7323ab0aad5SBarry Smith   PetscFunctionBegin;
7333ab0aad5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7343ab0aad5SBarry Smith   if (rhs) {
7353ab0aad5SBarry Smith     PetscValidHeaderSpecific(rhs,VEC_COOKIE,2);
7363ab0aad5SBarry Smith     PetscCheckSameComm(snes,1,rhs,2);
7373ab0aad5SBarry Smith     ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr);
7383ab0aad5SBarry Smith   }
7393ab0aad5SBarry Smith   if (snes->afine) {
7403ab0aad5SBarry Smith     ierr = VecDestroy(snes->afine);CHKERRQ(ierr);
7413ab0aad5SBarry Smith   }
7423ab0aad5SBarry Smith   snes->afine = rhs;
7433ab0aad5SBarry Smith   PetscFunctionReturn(0);
7443ab0aad5SBarry Smith }
7453ab0aad5SBarry Smith 
7464a2ae208SSatish Balay #undef __FUNCT__
7474a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
7489b94acceSBarry Smith /*@
74936851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
7509b94acceSBarry Smith                          SNESSetFunction().
7519b94acceSBarry Smith 
752c7afd0dbSLois Curfman McInnes    Collective on SNES
753c7afd0dbSLois Curfman McInnes 
7549b94acceSBarry Smith    Input Parameters:
755c7afd0dbSLois Curfman McInnes +  snes - the SNES context
756c7afd0dbSLois Curfman McInnes -  x - input vector
7579b94acceSBarry Smith 
7589b94acceSBarry Smith    Output Parameter:
7593638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
7609b94acceSBarry Smith 
7611bffabb2SLois Curfman McInnes    Notes:
76236851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
76336851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
76436851e7fSLois Curfman McInnes    themselves.
76536851e7fSLois Curfman McInnes 
76636851e7fSLois Curfman McInnes    Level: developer
76736851e7fSLois Curfman McInnes 
7689b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7699b94acceSBarry Smith 
770a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7719b94acceSBarry Smith @*/
7729b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x,Vec y)
7739b94acceSBarry Smith {
7749b94acceSBarry Smith   int    ierr;
7759b94acceSBarry Smith 
7763a40ed3dSBarry Smith   PetscFunctionBegin;
7774482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7784482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
7794482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
780c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
781c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,3);
782184914b5SBarry Smith 
783d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
784d64ed03dSBarry Smith   PetscStackPush("SNES user function");
7859b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
786d64ed03dSBarry Smith   PetscStackPop;
7873ab0aad5SBarry Smith   if (snes->afine) {
7883ab0aad5SBarry Smith     PetscScalar mone = -1.0;
7893ab0aad5SBarry Smith     ierr = VecAXPY(&mone,snes->afine,y);CHKERRQ(ierr);
7903ab0aad5SBarry Smith   }
791ae3c334cSLois Curfman McInnes   snes->nfuncs++;
792d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
7933a40ed3dSBarry Smith   PetscFunctionReturn(0);
7949b94acceSBarry Smith }
7959b94acceSBarry Smith 
7964a2ae208SSatish Balay #undef __FUNCT__
7974a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
79862fef451SLois Curfman McInnes /*@
79962fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
80062fef451SLois Curfman McInnes    set with SNESSetJacobian().
80162fef451SLois Curfman McInnes 
802c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
803c7afd0dbSLois Curfman McInnes 
80462fef451SLois Curfman McInnes    Input Parameters:
805c7afd0dbSLois Curfman McInnes +  snes - the SNES context
806c7afd0dbSLois Curfman McInnes -  x - input vector
80762fef451SLois Curfman McInnes 
80862fef451SLois Curfman McInnes    Output Parameters:
809c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
81062fef451SLois Curfman McInnes .  B - optional preconditioning matrix
811c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
812fee21e36SBarry Smith 
81362fef451SLois Curfman McInnes    Notes:
81462fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
81562fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
81662fef451SLois Curfman McInnes 
81794b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
818dc5a77f8SLois Curfman McInnes    flag parameter.
81962fef451SLois Curfman McInnes 
82036851e7fSLois Curfman McInnes    Level: developer
82136851e7fSLois Curfman McInnes 
82262fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
82362fef451SLois Curfman McInnes 
82494b7f48cSBarry Smith .seealso:  SNESSetJacobian(), KSPSetOperators()
82562fef451SLois Curfman McInnes @*/
8269b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8279b94acceSBarry Smith {
8289b94acceSBarry Smith   int    ierr;
8293a40ed3dSBarry Smith 
8303a40ed3dSBarry Smith   PetscFunctionBegin;
8314482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8324482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
8334482741eSBarry Smith   PetscValidPointer(flg,5);
834c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,X,2);
8353a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
836d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
837c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
838d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
8399b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
840d64ed03dSBarry Smith   PetscStackPop;
841d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
8426d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
8434482741eSBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
8444482741eSBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE,4);
8453a40ed3dSBarry Smith   PetscFunctionReturn(0);
8469b94acceSBarry Smith }
8479b94acceSBarry Smith 
8484a2ae208SSatish Balay #undef __FUNCT__
8494a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
8509b94acceSBarry Smith /*@C
8519b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
852044dda88SLois Curfman McInnes    location to store the matrix.
8539b94acceSBarry Smith 
854c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
855c7afd0dbSLois Curfman McInnes 
8569b94acceSBarry Smith    Input Parameters:
857c7afd0dbSLois Curfman McInnes +  snes - the SNES context
8589b94acceSBarry Smith .  A - Jacobian matrix
8599b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
8609b94acceSBarry Smith .  func - Jacobian evaluation routine
861c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
8622cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
8639b94acceSBarry Smith 
8649b94acceSBarry Smith    Calling sequence of func:
8658d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
8669b94acceSBarry Smith 
867c7afd0dbSLois Curfman McInnes +  x - input vector
8689b94acceSBarry Smith .  A - Jacobian matrix
8699b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
870ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
87194b7f48cSBarry Smith    structure (same as flag in KSPSetOperators())
872c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
8739b94acceSBarry Smith 
8749b94acceSBarry Smith    Notes:
87594b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
8762cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
877ac21db08SLois Curfman McInnes 
878ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
8799b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
8809b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
8819b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
8829b94acceSBarry Smith    throughout the global iterations.
8839b94acceSBarry Smith 
88436851e7fSLois Curfman McInnes    Level: beginner
88536851e7fSLois Curfman McInnes 
8869b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
8879b94acceSBarry Smith 
88894b7f48cSBarry Smith .seealso: KSPSetOperators(), SNESSetFunction()
8899b94acceSBarry Smith @*/
890454a90a3SBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
8919b94acceSBarry Smith {
8923a7fca6bSBarry Smith   int ierr;
8933a7fca6bSBarry Smith 
8943a40ed3dSBarry Smith   PetscFunctionBegin;
8954482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8964482741eSBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
8974482741eSBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
898c9780b6fSBarry Smith   if (A) PetscCheckSameComm(snes,1,A,2);
899c9780b6fSBarry Smith   if (B) PetscCheckSameComm(snes,1,B,2);
9003a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
9013a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
9023a7fca6bSBarry Smith   if (A) {
9033a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
9049b94acceSBarry Smith     snes->jacobian = A;
9053a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
9063a7fca6bSBarry Smith   }
9073a7fca6bSBarry Smith   if (B) {
9083a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
9099b94acceSBarry Smith     snes->jacobian_pre = B;
9103a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
9113a7fca6bSBarry Smith   }
9123a40ed3dSBarry Smith   PetscFunctionReturn(0);
9139b94acceSBarry Smith }
91462fef451SLois Curfman McInnes 
9154a2ae208SSatish Balay #undef __FUNCT__
9164a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
917c2aafc4cSSatish Balay /*@C
918b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
919b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
920b4fd4287SBarry Smith 
921c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
922c7afd0dbSLois Curfman McInnes 
923b4fd4287SBarry Smith    Input Parameter:
924b4fd4287SBarry Smith .  snes - the nonlinear solver context
925b4fd4287SBarry Smith 
926b4fd4287SBarry Smith    Output Parameters:
927c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
928b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
92900036973SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
93000036973SBarry Smith -  func - location to put Jacobian function (or PETSC_NULL)
931fee21e36SBarry Smith 
93236851e7fSLois Curfman McInnes    Level: advanced
93336851e7fSLois Curfman McInnes 
934b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
935b4fd4287SBarry Smith @*/
93600036973SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
937b4fd4287SBarry Smith {
9383a40ed3dSBarry Smith   PetscFunctionBegin;
9394482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
940b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
941b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
942b4fd4287SBarry Smith   if (ctx)  *ctx  = snes->jacP;
94300036973SBarry Smith   if (func) *func = snes->computejacobian;
9443a40ed3dSBarry Smith   PetscFunctionReturn(0);
945b4fd4287SBarry Smith }
946b4fd4287SBarry Smith 
9479b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
94845fc7adcSBarry Smith extern int SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
9499b94acceSBarry Smith 
9504a2ae208SSatish Balay #undef __FUNCT__
9514a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
9529b94acceSBarry Smith /*@
9539b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
954272ac6f2SLois Curfman McInnes    of a nonlinear solver.
9559b94acceSBarry Smith 
956fee21e36SBarry Smith    Collective on SNES
957fee21e36SBarry Smith 
958c7afd0dbSLois Curfman McInnes    Input Parameters:
959c7afd0dbSLois Curfman McInnes +  snes - the SNES context
960c7afd0dbSLois Curfman McInnes -  x - the solution vector
961c7afd0dbSLois Curfman McInnes 
962272ac6f2SLois Curfman McInnes    Notes:
963272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
964272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
965272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
966272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
967272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
968272ac6f2SLois Curfman McInnes 
96936851e7fSLois Curfman McInnes    Level: advanced
97036851e7fSLois Curfman McInnes 
9719b94acceSBarry Smith .keywords: SNES, nonlinear, setup
9729b94acceSBarry Smith 
9739b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
9749b94acceSBarry Smith @*/
9758ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x)
9769b94acceSBarry Smith {
977f1af5d2fSBarry Smith   int        ierr;
9784b27c08aSLois Curfman McInnes   PetscTruth flg, iseqtr;
9793a40ed3dSBarry Smith 
9803a40ed3dSBarry Smith   PetscFunctionBegin;
9814482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9824482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
983c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
9848ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
9859b94acceSBarry Smith 
986b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
987c1f60f51SBarry Smith   /*
988c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
989dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
990c1f60f51SBarry Smith   */
991c1f60f51SBarry Smith   if (flg) {
992c1f60f51SBarry Smith     Mat J;
9935a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
9945a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
9953a7fca6bSBarry Smith     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n");
9963a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
9973a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
998c1f60f51SBarry Smith   }
99945fc7adcSBarry Smith 
1000b4014bb2SMatthew Knepley #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
100145fc7adcSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
100245fc7adcSBarry Smith   if (flg) {
100345fc7adcSBarry Smith     Mat J;
100445fc7adcSBarry Smith     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
100545fc7adcSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
100645fc7adcSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
100745fc7adcSBarry Smith   }
100832a4b47aSMatthew Knepley #endif
100945fc7adcSBarry Smith 
1010b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1011c1f60f51SBarry Smith   /*
1012dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1013c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1014c1f60f51SBarry Smith    */
101531615d04SLois Curfman McInnes   if (flg) {
1016272ac6f2SLois Curfman McInnes     Mat  J;
1017b5d62d44SBarry Smith     KSP ksp;
101894b7f48cSBarry Smith     PC   pc;
1019f3ef73ceSBarry Smith 
10205a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10213a7fca6bSBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
102293b98531SBarry Smith     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n");
10233a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
10243a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
10253a7fca6bSBarry Smith 
1026f3ef73ceSBarry Smith     /* force no preconditioner */
102794b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1028b5d62d44SBarry Smith     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1029a9815358SBarry Smith     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1030a9815358SBarry Smith     if (!flg) {
1031f3ef73ceSBarry Smith       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1032272ac6f2SLois Curfman McInnes     }
1033a9815358SBarry Smith   }
1034f3ef73ceSBarry Smith 
103529bbc08cSBarry Smith   if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
103629bbc08cSBarry Smith   if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
103729bbc08cSBarry Smith   if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1038a8c6a408SBarry Smith   if (snes->vec_func == snes->vec_sol) {
103929bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1040a8c6a408SBarry Smith   }
1041a703fe33SLois Curfman McInnes 
1042a703fe33SLois Curfman McInnes   /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
10434b27c08aSLois Curfman McInnes   ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr);
10446831982aSBarry Smith   if (snes->ksp_ewconv && !iseqtr) {
104594b7f48cSBarry Smith     KSP ksp;
104694b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1047186905e3SBarry Smith     ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1048a703fe33SLois Curfman McInnes   }
10494b27c08aSLois Curfman McInnes 
1050a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
105182bf6240SBarry Smith   snes->setupcalled = 1;
10523a40ed3dSBarry Smith   PetscFunctionReturn(0);
10539b94acceSBarry Smith }
10549b94acceSBarry Smith 
10554a2ae208SSatish Balay #undef __FUNCT__
10564a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
10579b94acceSBarry Smith /*@C
10589b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
10599b94acceSBarry Smith    with SNESCreate().
10609b94acceSBarry Smith 
1061c7afd0dbSLois Curfman McInnes    Collective on SNES
1062c7afd0dbSLois Curfman McInnes 
10639b94acceSBarry Smith    Input Parameter:
10649b94acceSBarry Smith .  snes - the SNES context
10659b94acceSBarry Smith 
106636851e7fSLois Curfman McInnes    Level: beginner
106736851e7fSLois Curfman McInnes 
10689b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
10699b94acceSBarry Smith 
107063a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
10719b94acceSBarry Smith @*/
10729b94acceSBarry Smith int SNESDestroy(SNES snes)
10739b94acceSBarry Smith {
1074b8a78c4aSBarry Smith   int i,ierr;
10753a40ed3dSBarry Smith 
10763a40ed3dSBarry Smith   PetscFunctionBegin;
10774482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
10783a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1079d4bb536fSBarry Smith 
1080be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
10810f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1082be0abb6dSBarry Smith 
1083e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1084606d414cSSatish Balay   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
10853a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
10863a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
10873ab0aad5SBarry Smith   if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);}
108894b7f48cSBarry Smith   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
1089522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1090b8a78c4aSBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1091b8a78c4aSBarry Smith     if (snes->monitordestroy[i]) {
1092b8a78c4aSBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1093b8a78c4aSBarry Smith     }
1094b8a78c4aSBarry Smith   }
1095b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)snes);
10960452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
10973a40ed3dSBarry Smith   PetscFunctionReturn(0);
10989b94acceSBarry Smith }
10999b94acceSBarry Smith 
11009b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11019b94acceSBarry Smith 
11024a2ae208SSatish Balay #undef __FUNCT__
11034a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
11049b94acceSBarry Smith /*@
1105d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11069b94acceSBarry Smith 
1107c7afd0dbSLois Curfman McInnes    Collective on SNES
1108c7afd0dbSLois Curfman McInnes 
11099b94acceSBarry Smith    Input Parameters:
1110c7afd0dbSLois Curfman McInnes +  snes - the SNES context
111133174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
111233174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
111333174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
111433174efeSLois Curfman McInnes            of the change in the solution between steps
111533174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1116c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1117fee21e36SBarry Smith 
111833174efeSLois Curfman McInnes    Options Database Keys:
1119c7afd0dbSLois Curfman McInnes +    -snes_atol <atol> - Sets atol
1120c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1121c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1122c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1123c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
11249b94acceSBarry Smith 
1125d7a720efSLois Curfman McInnes    Notes:
11269b94acceSBarry Smith    The default maximum number of iterations is 50.
11279b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
11289b94acceSBarry Smith 
112936851e7fSLois Curfman McInnes    Level: intermediate
113036851e7fSLois Curfman McInnes 
113133174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
11329b94acceSBarry Smith 
1133d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
11349b94acceSBarry Smith @*/
1135329f5518SBarry Smith int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
11369b94acceSBarry Smith {
11373a40ed3dSBarry Smith   PetscFunctionBegin;
11384482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1139d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1140d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1141d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1142d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1143d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
11443a40ed3dSBarry Smith   PetscFunctionReturn(0);
11459b94acceSBarry Smith }
11469b94acceSBarry Smith 
11474a2ae208SSatish Balay #undef __FUNCT__
11484a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
11499b94acceSBarry Smith /*@
115033174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
115133174efeSLois Curfman McInnes 
1152c7afd0dbSLois Curfman McInnes    Not Collective
1153c7afd0dbSLois Curfman McInnes 
115433174efeSLois Curfman McInnes    Input Parameters:
1155c7afd0dbSLois Curfman McInnes +  snes - the SNES context
115633174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
115733174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
115833174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
115933174efeSLois Curfman McInnes            of the change in the solution between steps
116033174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1161c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1162fee21e36SBarry Smith 
116333174efeSLois Curfman McInnes    Notes:
116433174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
116533174efeSLois Curfman McInnes 
116636851e7fSLois Curfman McInnes    Level: intermediate
116736851e7fSLois Curfman McInnes 
116833174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
116933174efeSLois Curfman McInnes 
117033174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
117133174efeSLois Curfman McInnes @*/
1172329f5518SBarry Smith int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
117333174efeSLois Curfman McInnes {
11743a40ed3dSBarry Smith   PetscFunctionBegin;
11754482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
117633174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
117733174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
117833174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
117933174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
118033174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
11813a40ed3dSBarry Smith   PetscFunctionReturn(0);
118233174efeSLois Curfman McInnes }
118333174efeSLois Curfman McInnes 
11844a2ae208SSatish Balay #undef __FUNCT__
11854a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
118633174efeSLois Curfman McInnes /*@
11879b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
11889b94acceSBarry Smith 
1189fee21e36SBarry Smith    Collective on SNES
1190fee21e36SBarry Smith 
1191c7afd0dbSLois Curfman McInnes    Input Parameters:
1192c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1193c7afd0dbSLois Curfman McInnes -  tol - tolerance
1194c7afd0dbSLois Curfman McInnes 
11959b94acceSBarry Smith    Options Database Key:
1196c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
11979b94acceSBarry Smith 
119836851e7fSLois Curfman McInnes    Level: intermediate
119936851e7fSLois Curfman McInnes 
12009b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12019b94acceSBarry Smith 
1202d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
12039b94acceSBarry Smith @*/
1204329f5518SBarry Smith int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
12059b94acceSBarry Smith {
12063a40ed3dSBarry Smith   PetscFunctionBegin;
12074482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
12089b94acceSBarry Smith   snes->deltatol = tol;
12093a40ed3dSBarry Smith   PetscFunctionReturn(0);
12109b94acceSBarry Smith }
12119b94acceSBarry Smith 
1212df9fa365SBarry Smith /*
1213df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1214df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1215df9fa365SBarry Smith    macros instead of functions
1216df9fa365SBarry Smith */
12174a2ae208SSatish Balay #undef __FUNCT__
12184a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor"
1219329f5518SBarry Smith int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1220ce1608b8SBarry Smith {
1221ce1608b8SBarry Smith   int ierr;
1222ce1608b8SBarry Smith 
1223ce1608b8SBarry Smith   PetscFunctionBegin;
12244482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1225ce1608b8SBarry Smith   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1226ce1608b8SBarry Smith   PetscFunctionReturn(0);
1227ce1608b8SBarry Smith }
1228ce1608b8SBarry Smith 
12294a2ae208SSatish Balay #undef __FUNCT__
12304a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate"
12311836bdbcSSatish Balay int SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1232df9fa365SBarry Smith {
1233df9fa365SBarry Smith   int ierr;
1234df9fa365SBarry Smith 
1235df9fa365SBarry Smith   PetscFunctionBegin;
1236df9fa365SBarry Smith   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1237df9fa365SBarry Smith   PetscFunctionReturn(0);
1238df9fa365SBarry Smith }
1239df9fa365SBarry Smith 
12404a2ae208SSatish Balay #undef __FUNCT__
12414a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy"
1242b0a32e0cSBarry Smith int SNESLGMonitorDestroy(PetscDrawLG draw)
1243df9fa365SBarry Smith {
1244df9fa365SBarry Smith   int ierr;
1245df9fa365SBarry Smith 
1246df9fa365SBarry Smith   PetscFunctionBegin;
1247df9fa365SBarry Smith   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1248df9fa365SBarry Smith   PetscFunctionReturn(0);
1249df9fa365SBarry Smith }
1250df9fa365SBarry Smith 
12519b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
12529b94acceSBarry Smith 
12534a2ae208SSatish Balay #undef __FUNCT__
12544a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor"
12559b94acceSBarry Smith /*@C
12565cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
12579b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
12589b94acceSBarry Smith    progress.
12599b94acceSBarry Smith 
1260fee21e36SBarry Smith    Collective on SNES
1261fee21e36SBarry Smith 
1262c7afd0dbSLois Curfman McInnes    Input Parameters:
1263c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1264c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1265b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1266b3006f0bSLois Curfman McInnes           monitor routine (use PETSC_NULL if no context is desitre)
1267b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1268b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
12699b94acceSBarry Smith 
1270c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1271329f5518SBarry Smith $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1272c7afd0dbSLois Curfman McInnes 
1273c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1274c7afd0dbSLois Curfman McInnes .    its - iteration number
1275c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
127640a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
12779b94acceSBarry Smith 
12789665c990SLois Curfman McInnes    Options Database Keys:
1279c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1280c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1281c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1282c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1283c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1284c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1285c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1286c7afd0dbSLois Curfman McInnes                             the options database.
12879665c990SLois Curfman McInnes 
1288639f9d9dSBarry Smith    Notes:
12896bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
12906bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
12916bc08f3fSLois Curfman McInnes    order in which they were set.
1292639f9d9dSBarry Smith 
129336851e7fSLois Curfman McInnes    Level: intermediate
129436851e7fSLois Curfman McInnes 
12959b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
12969b94acceSBarry Smith 
12975cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
12989b94acceSBarry Smith @*/
1299329f5518SBarry Smith int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *))
13009b94acceSBarry Smith {
13013a40ed3dSBarry Smith   PetscFunctionBegin;
13024482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1303639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
130429bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1305639f9d9dSBarry Smith   }
1306639f9d9dSBarry Smith 
1307639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1308b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1309639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
13103a40ed3dSBarry Smith   PetscFunctionReturn(0);
13119b94acceSBarry Smith }
13129b94acceSBarry Smith 
13134a2ae208SSatish Balay #undef __FUNCT__
13144a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor"
13155cd90555SBarry Smith /*@C
13165cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
13175cd90555SBarry Smith 
1318c7afd0dbSLois Curfman McInnes    Collective on SNES
1319c7afd0dbSLois Curfman McInnes 
13205cd90555SBarry Smith    Input Parameters:
13215cd90555SBarry Smith .  snes - the SNES context
13225cd90555SBarry Smith 
13231a480d89SAdministrator    Options Database Key:
1324c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1325c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1326c7afd0dbSLois Curfman McInnes     set via the options database
13275cd90555SBarry Smith 
13285cd90555SBarry Smith    Notes:
13295cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
13305cd90555SBarry Smith 
133136851e7fSLois Curfman McInnes    Level: intermediate
133236851e7fSLois Curfman McInnes 
13335cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
13345cd90555SBarry Smith 
13355cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
13365cd90555SBarry Smith @*/
13375cd90555SBarry Smith int SNESClearMonitor(SNES snes)
13385cd90555SBarry Smith {
13395cd90555SBarry Smith   PetscFunctionBegin;
13404482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
13415cd90555SBarry Smith   snes->numbermonitors = 0;
13425cd90555SBarry Smith   PetscFunctionReturn(0);
13435cd90555SBarry Smith }
13445cd90555SBarry Smith 
13454a2ae208SSatish Balay #undef __FUNCT__
13464a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
13479b94acceSBarry Smith /*@C
13489b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
13499b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
13509b94acceSBarry Smith 
1351fee21e36SBarry Smith    Collective on SNES
1352fee21e36SBarry Smith 
1353c7afd0dbSLois Curfman McInnes    Input Parameters:
1354c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1355c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1356c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1357c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
13589b94acceSBarry Smith 
1359c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1360329f5518SBarry Smith $     int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1361c7afd0dbSLois Curfman McInnes 
1362c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1363c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1364184914b5SBarry Smith .    reason - reason for convergence/divergence
1365c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
13664b27c08aSLois Curfman McInnes .    gnorm - 2-norm of current step
13674b27c08aSLois Curfman McInnes -    f - 2-norm of function
13689b94acceSBarry Smith 
136936851e7fSLois Curfman McInnes    Level: advanced
137036851e7fSLois Curfman McInnes 
13719b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
13729b94acceSBarry Smith 
13734b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR()
13749b94acceSBarry Smith @*/
1375329f5518SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
13769b94acceSBarry Smith {
13773a40ed3dSBarry Smith   PetscFunctionBegin;
13784482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
13799b94acceSBarry Smith   (snes)->converged = func;
13809b94acceSBarry Smith   (snes)->cnvP      = cctx;
13813a40ed3dSBarry Smith   PetscFunctionReturn(0);
13829b94acceSBarry Smith }
13839b94acceSBarry Smith 
13844a2ae208SSatish Balay #undef __FUNCT__
13854a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
1386184914b5SBarry Smith /*@C
1387184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1388184914b5SBarry Smith 
1389184914b5SBarry Smith    Not Collective
1390184914b5SBarry Smith 
1391184914b5SBarry Smith    Input Parameter:
1392184914b5SBarry Smith .  snes - the SNES context
1393184914b5SBarry Smith 
1394184914b5SBarry Smith    Output Parameter:
1395e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1396184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1397184914b5SBarry Smith 
1398184914b5SBarry Smith    Level: intermediate
1399184914b5SBarry Smith 
1400184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1401184914b5SBarry Smith 
1402184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1403184914b5SBarry Smith 
14044b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason
1405184914b5SBarry Smith @*/
1406184914b5SBarry Smith int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1407184914b5SBarry Smith {
1408184914b5SBarry Smith   PetscFunctionBegin;
14094482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14104482741eSBarry Smith   PetscValidPointer(reason,2);
1411184914b5SBarry Smith   *reason = snes->reason;
1412184914b5SBarry Smith   PetscFunctionReturn(0);
1413184914b5SBarry Smith }
1414184914b5SBarry Smith 
14154a2ae208SSatish Balay #undef __FUNCT__
14164a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1417c9005455SLois Curfman McInnes /*@
1418c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1419c9005455SLois Curfman McInnes 
1420fee21e36SBarry Smith    Collective on SNES
1421fee21e36SBarry Smith 
1422c7afd0dbSLois Curfman McInnes    Input Parameters:
1423c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1424c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1425cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1426758f92a0SBarry Smith .  na  - size of a and its
142764731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1428758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1429c7afd0dbSLois Curfman McInnes 
1430c9005455SLois Curfman McInnes    Notes:
14314b27c08aSLois Curfman McInnes    If set, this array will contain the function norms computed
1432c9005455SLois Curfman McInnes    at each step.
1433c9005455SLois Curfman McInnes 
1434c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1435c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1436c9005455SLois Curfman McInnes    during the section of code that is being timed.
1437c9005455SLois Curfman McInnes 
143836851e7fSLois Curfman McInnes    Level: intermediate
143936851e7fSLois Curfman McInnes 
1440c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1441758f92a0SBarry Smith 
144208405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1443758f92a0SBarry Smith 
1444c9005455SLois Curfman McInnes @*/
14451836bdbcSSatish Balay int SNESSetConvergenceHistory(SNES snes,PetscReal a[],int *its,int na,PetscTruth reset)
1446c9005455SLois Curfman McInnes {
14473a40ed3dSBarry Smith   PetscFunctionBegin;
14484482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14494482741eSBarry Smith   if (na) PetscValidScalarPointer(a,2);
1450c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1451758f92a0SBarry Smith   snes->conv_hist_its   = its;
1452758f92a0SBarry Smith   snes->conv_hist_max   = na;
1453758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1454758f92a0SBarry Smith   PetscFunctionReturn(0);
1455758f92a0SBarry Smith }
1456758f92a0SBarry Smith 
14574a2ae208SSatish Balay #undef __FUNCT__
14584a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
14590c4c9dddSBarry Smith /*@C
1460758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1461758f92a0SBarry Smith 
1462758f92a0SBarry Smith    Collective on SNES
1463758f92a0SBarry Smith 
1464758f92a0SBarry Smith    Input Parameter:
1465758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1466758f92a0SBarry Smith 
1467758f92a0SBarry Smith    Output Parameters:
1468758f92a0SBarry Smith .  a   - array to hold history
1469758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1470758f92a0SBarry Smith          negative if not converged) for each solve.
1471758f92a0SBarry Smith -  na  - size of a and its
1472758f92a0SBarry Smith 
1473758f92a0SBarry Smith    Notes:
1474758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1475758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1476758f92a0SBarry Smith 
1477758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1478758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1479758f92a0SBarry Smith    during the section of code that is being timed.
1480758f92a0SBarry Smith 
1481758f92a0SBarry Smith    Level: intermediate
1482758f92a0SBarry Smith 
1483758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1484758f92a0SBarry Smith 
1485758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1486758f92a0SBarry Smith 
1487758f92a0SBarry Smith @*/
14881836bdbcSSatish Balay int SNESGetConvergenceHistory(SNES snes,PetscReal *a[],int *its[],int *na)
1489758f92a0SBarry Smith {
1490758f92a0SBarry Smith   PetscFunctionBegin;
14914482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1492758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1493758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1494758f92a0SBarry Smith   if (na) *na   = snes->conv_hist_len;
14953a40ed3dSBarry Smith   PetscFunctionReturn(0);
1496c9005455SLois Curfman McInnes }
1497c9005455SLois Curfman McInnes 
1498e74ef692SMatthew Knepley #undef __FUNCT__
1499e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC"
1500ac226902SBarry Smith /*@C
150176b2cf59SMatthew Knepley   SNESSetRhsBC - Sets the function which applies boundary conditions
150276b2cf59SMatthew Knepley   to the Rhs of each system.
150376b2cf59SMatthew Knepley 
150476b2cf59SMatthew Knepley   Collective on SNES
150576b2cf59SMatthew Knepley 
150676b2cf59SMatthew Knepley   Input Parameters:
150776b2cf59SMatthew Knepley . snes - The nonlinear solver context
150876b2cf59SMatthew Knepley . func - The function
150976b2cf59SMatthew Knepley 
151076b2cf59SMatthew Knepley   Calling sequence of func:
151176b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx);
151276b2cf59SMatthew Knepley 
151376b2cf59SMatthew Knepley . rhs - The current rhs vector
151476b2cf59SMatthew Knepley . ctx - The user-context
151576b2cf59SMatthew Knepley 
151676b2cf59SMatthew Knepley   Level: intermediate
151776b2cf59SMatthew Knepley 
151876b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions
151976b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate()
152076b2cf59SMatthew Knepley @*/
152176b2cf59SMatthew Knepley int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *))
152276b2cf59SMatthew Knepley {
152376b2cf59SMatthew Knepley   PetscFunctionBegin;
15244482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
152576b2cf59SMatthew Knepley   snes->applyrhsbc = func;
152676b2cf59SMatthew Knepley   PetscFunctionReturn(0);
152776b2cf59SMatthew Knepley }
152876b2cf59SMatthew Knepley 
1529e74ef692SMatthew Knepley #undef __FUNCT__
1530e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC"
153176b2cf59SMatthew Knepley /*@
153276b2cf59SMatthew Knepley   SNESDefaultRhsBC - The default boundary condition function which does nothing.
153376b2cf59SMatthew Knepley 
153476b2cf59SMatthew Knepley   Not collective
153576b2cf59SMatthew Knepley 
153676b2cf59SMatthew Knepley   Input Parameters:
153776b2cf59SMatthew Knepley . snes - The nonlinear solver context
153876b2cf59SMatthew Knepley . rhs  - The Rhs
153976b2cf59SMatthew Knepley . ctx  - The user-context
154076b2cf59SMatthew Knepley 
154176b2cf59SMatthew Knepley   Level: beginner
154276b2cf59SMatthew Knepley 
154376b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions
154476b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate()
154576b2cf59SMatthew Knepley @*/
154676b2cf59SMatthew Knepley int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx)
154776b2cf59SMatthew Knepley {
154876b2cf59SMatthew Knepley   PetscFunctionBegin;
154976b2cf59SMatthew Knepley   PetscFunctionReturn(0);
155076b2cf59SMatthew Knepley }
155176b2cf59SMatthew Knepley 
1552e74ef692SMatthew Knepley #undef __FUNCT__
1553e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC"
1554ac226902SBarry Smith /*@C
155576b2cf59SMatthew Knepley   SNESSetSolutionBC - Sets the function which applies boundary conditions
155676b2cf59SMatthew Knepley   to the solution of each system.
155776b2cf59SMatthew Knepley 
155876b2cf59SMatthew Knepley   Collective on SNES
155976b2cf59SMatthew Knepley 
156076b2cf59SMatthew Knepley   Input Parameters:
156176b2cf59SMatthew Knepley . snes - The nonlinear solver context
156276b2cf59SMatthew Knepley . func - The function
156376b2cf59SMatthew Knepley 
156476b2cf59SMatthew Knepley   Calling sequence of func:
15659d06fb6bSMatthew Knepley . func (SNES snes, Vec rsol, void *ctx);
156676b2cf59SMatthew Knepley 
156776b2cf59SMatthew Knepley . sol - The current solution vector
156876b2cf59SMatthew Knepley . ctx - The user-context
156976b2cf59SMatthew Knepley 
157076b2cf59SMatthew Knepley   Level: intermediate
157176b2cf59SMatthew Knepley 
157276b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions
157376b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate()
157476b2cf59SMatthew Knepley @*/
157576b2cf59SMatthew Knepley int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *))
157676b2cf59SMatthew Knepley {
157776b2cf59SMatthew Knepley   PetscFunctionBegin;
15784482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
157976b2cf59SMatthew Knepley   snes->applysolbc = func;
158076b2cf59SMatthew Knepley   PetscFunctionReturn(0);
158176b2cf59SMatthew Knepley }
158276b2cf59SMatthew Knepley 
1583e74ef692SMatthew Knepley #undef __FUNCT__
1584e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC"
158576b2cf59SMatthew Knepley /*@
158676b2cf59SMatthew Knepley   SNESDefaultSolutionBC - The default boundary condition function which does nothing.
158776b2cf59SMatthew Knepley 
158876b2cf59SMatthew Knepley   Not collective
158976b2cf59SMatthew Knepley 
159076b2cf59SMatthew Knepley   Input Parameters:
159176b2cf59SMatthew Knepley . snes - The nonlinear solver context
159276b2cf59SMatthew Knepley . sol  - The solution
159376b2cf59SMatthew Knepley . ctx  - The user-context
159476b2cf59SMatthew Knepley 
159576b2cf59SMatthew Knepley   Level: beginner
159676b2cf59SMatthew Knepley 
159776b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions
159876b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate()
159976b2cf59SMatthew Knepley @*/
160076b2cf59SMatthew Knepley int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx)
160176b2cf59SMatthew Knepley {
160276b2cf59SMatthew Knepley   PetscFunctionBegin;
160376b2cf59SMatthew Knepley   PetscFunctionReturn(0);
160476b2cf59SMatthew Knepley }
160576b2cf59SMatthew Knepley 
1606e74ef692SMatthew Knepley #undef __FUNCT__
1607e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
1608ac226902SBarry Smith /*@C
160976b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
161076b2cf59SMatthew Knepley   at the beginning of every step of the iteration.
161176b2cf59SMatthew Knepley 
161276b2cf59SMatthew Knepley   Collective on SNES
161376b2cf59SMatthew Knepley 
161476b2cf59SMatthew Knepley   Input Parameters:
161576b2cf59SMatthew Knepley . snes - The nonlinear solver context
161676b2cf59SMatthew Knepley . func - The function
161776b2cf59SMatthew Knepley 
161876b2cf59SMatthew Knepley   Calling sequence of func:
161976b2cf59SMatthew Knepley . func (TS ts, int step);
162076b2cf59SMatthew Knepley 
162176b2cf59SMatthew Knepley . step - The current step of the iteration
162276b2cf59SMatthew Knepley 
162376b2cf59SMatthew Knepley   Level: intermediate
162476b2cf59SMatthew Knepley 
162576b2cf59SMatthew Knepley .keywords: SNES, update
162676b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
162776b2cf59SMatthew Knepley @*/
162876b2cf59SMatthew Knepley int SNESSetUpdate(SNES snes, int (*func)(SNES, int))
162976b2cf59SMatthew Knepley {
163076b2cf59SMatthew Knepley   PetscFunctionBegin;
16314482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
163276b2cf59SMatthew Knepley   snes->update = func;
163376b2cf59SMatthew Knepley   PetscFunctionReturn(0);
163476b2cf59SMatthew Knepley }
163576b2cf59SMatthew Knepley 
1636e74ef692SMatthew Knepley #undef __FUNCT__
1637e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
163876b2cf59SMatthew Knepley /*@
163976b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
164076b2cf59SMatthew Knepley 
164176b2cf59SMatthew Knepley   Not collective
164276b2cf59SMatthew Knepley 
164376b2cf59SMatthew Knepley   Input Parameters:
164476b2cf59SMatthew Knepley . snes - The nonlinear solver context
164576b2cf59SMatthew Knepley . step - The current step of the iteration
164676b2cf59SMatthew Knepley 
1647205452f4SMatthew Knepley   Level: intermediate
1648205452f4SMatthew Knepley 
164976b2cf59SMatthew Knepley .keywords: SNES, update
165076b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
165176b2cf59SMatthew Knepley @*/
165276b2cf59SMatthew Knepley int SNESDefaultUpdate(SNES snes, int step)
165376b2cf59SMatthew Knepley {
165476b2cf59SMatthew Knepley   PetscFunctionBegin;
165576b2cf59SMatthew Knepley   PetscFunctionReturn(0);
165676b2cf59SMatthew Knepley }
165776b2cf59SMatthew Knepley 
16584a2ae208SSatish Balay #undef __FUNCT__
16594a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
16609b94acceSBarry Smith /*
16619b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
16629b94acceSBarry Smith    positive parameter delta.
16639b94acceSBarry Smith 
16649b94acceSBarry Smith     Input Parameters:
1665c7afd0dbSLois Curfman McInnes +   snes - the SNES context
16669b94acceSBarry Smith .   y - approximate solution of linear system
16679b94acceSBarry Smith .   fnorm - 2-norm of current function
1668c7afd0dbSLois Curfman McInnes -   delta - trust region size
16699b94acceSBarry Smith 
16709b94acceSBarry Smith     Output Parameters:
1671c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
16729b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
16739b94acceSBarry Smith     region, and exceeds zero otherwise.
1674c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
16759b94acceSBarry Smith 
16769b94acceSBarry Smith     Note:
16774b27c08aSLois Curfman McInnes     For non-trust region methods such as SNESLS, the parameter delta
16789b94acceSBarry Smith     is set to be the maximum allowable step size.
16799b94acceSBarry Smith 
16809b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
16819b94acceSBarry Smith */
1682064f8208SBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
16839b94acceSBarry Smith {
1684064f8208SBarry Smith   PetscReal   nrm;
1685ea709b57SSatish Balay   PetscScalar cnorm;
16863a40ed3dSBarry Smith   int         ierr;
16873a40ed3dSBarry Smith 
16883a40ed3dSBarry Smith   PetscFunctionBegin;
16894482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16904482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1691c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,2);
1692184914b5SBarry Smith 
1693064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1694064f8208SBarry Smith   if (nrm > *delta) {
1695064f8208SBarry Smith      nrm = *delta/nrm;
1696064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
1697064f8208SBarry Smith      cnorm = nrm;
1698329f5518SBarry Smith      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
16999b94acceSBarry Smith      *ynorm = *delta;
17009b94acceSBarry Smith   } else {
17019b94acceSBarry Smith      *gpnorm = 0.0;
1702064f8208SBarry Smith      *ynorm = nrm;
17039b94acceSBarry Smith   }
17043a40ed3dSBarry Smith   PetscFunctionReturn(0);
17059b94acceSBarry Smith }
17069b94acceSBarry Smith 
17075968eb51SBarry Smith static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero",
17085968eb51SBarry Smith                                          "not currently used",
17095968eb51SBarry Smith                                          "line search failed",
17105968eb51SBarry Smith                                          "reach maximum number of iterations",
17115968eb51SBarry Smith                                          "function norm became NaN (not a number)",
17125968eb51SBarry Smith                                          "not currently used",
17135968eb51SBarry Smith                                          "number of function computations exceeded",
17145968eb51SBarry Smith                                          "not currently used",
17155968eb51SBarry Smith                                          "still iterating",
17165968eb51SBarry Smith                                          "not currently used",
17175968eb51SBarry Smith                                          "absolute size of function norm",
17185968eb51SBarry Smith                                          "relative decrease in function norm",
17195968eb51SBarry Smith                                          "step size is small",
17205968eb51SBarry Smith                                          "not currently used",
17215968eb51SBarry Smith                                          "not currently used",
17225968eb51SBarry Smith                                          "small trust region"};
17235968eb51SBarry Smith 
17244a2ae208SSatish Balay #undef __FUNCT__
17254a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
17269b94acceSBarry Smith /*@
17279b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
172863a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
17299b94acceSBarry Smith 
1730c7afd0dbSLois Curfman McInnes    Collective on SNES
1731c7afd0dbSLois Curfman McInnes 
1732b2002411SLois Curfman McInnes    Input Parameters:
1733c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1734c7afd0dbSLois Curfman McInnes -  x - the solution vector
17359b94acceSBarry Smith 
1736b2002411SLois Curfman McInnes    Notes:
17378ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
17388ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
17398ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
17408ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
17418ddd3da0SLois Curfman McInnes 
174236851e7fSLois Curfman McInnes    Level: beginner
174336851e7fSLois Curfman McInnes 
17449b94acceSBarry Smith .keywords: SNES, nonlinear, solve
17459b94acceSBarry Smith 
17463ab0aad5SBarry Smith .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs()
17479b94acceSBarry Smith @*/
1748c9780b6fSBarry Smith int SNESSolve(SNES snes,Vec x)
17499b94acceSBarry Smith {
1750f1af5d2fSBarry Smith   int        ierr;
1751f1af5d2fSBarry Smith   PetscTruth flg;
1752052efed2SBarry Smith 
17533a40ed3dSBarry Smith   PetscFunctionBegin;
17544482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
17554482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1756c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
175729bbc08cSBarry Smith   if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
175874637425SBarry Smith 
175982bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
1760c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
1761758f92a0SBarry Smith   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
1762d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
176350ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1764c9780b6fSBarry Smith   ierr = (*(snes)->solve)(snes);CHKERRQ(ierr);
1765d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1766b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
17678b179ff8SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); }
1768da9b6338SBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
1769da9b6338SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
17705968eb51SBarry Smith   if (snes->printreason) {
17715968eb51SBarry Smith     if (snes->reason > 0) {
17725968eb51SBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr);
17735968eb51SBarry Smith     } else {
17745968eb51SBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr);
17755968eb51SBarry Smith     }
17765968eb51SBarry Smith   }
17775968eb51SBarry Smith 
17783a40ed3dSBarry Smith   PetscFunctionReturn(0);
17799b94acceSBarry Smith }
17809b94acceSBarry Smith 
17819b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
17829b94acceSBarry Smith 
17834a2ae208SSatish Balay #undef __FUNCT__
17844a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
178582bf6240SBarry Smith /*@C
17864b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
17879b94acceSBarry Smith 
1788fee21e36SBarry Smith    Collective on SNES
1789fee21e36SBarry Smith 
1790c7afd0dbSLois Curfman McInnes    Input Parameters:
1791c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1792454a90a3SBarry Smith -  type - a known method
1793c7afd0dbSLois Curfman McInnes 
1794c7afd0dbSLois Curfman McInnes    Options Database Key:
1795454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
1796c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1797ae12b187SLois Curfman McInnes 
17989b94acceSBarry Smith    Notes:
1799e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
18004b27c08aSLois Curfman McInnes +    SNESLS - Newton's method with line search
1801c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
18024b27c08aSLois Curfman McInnes .    SNESTR - Newton's method with trust region
1803c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
18049b94acceSBarry Smith 
1805ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1806ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1807ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1808ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1809ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1810ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1811ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1812ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1813ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1814b0a32e0cSBarry Smith   appropriate method.
181536851e7fSLois Curfman McInnes 
181636851e7fSLois Curfman McInnes   Level: intermediate
1817a703fe33SLois Curfman McInnes 
1818454a90a3SBarry Smith .keywords: SNES, set, type
1819435da068SBarry Smith 
1820435da068SBarry Smith .seealso: SNESType, SNESCreate()
1821435da068SBarry Smith 
18229b94acceSBarry Smith @*/
18230e33f6ddSBarry Smith int SNESSetType(SNES snes,const SNESType type)
18249b94acceSBarry Smith {
18250f5bd95cSBarry Smith   int        ierr,(*r)(SNES);
18266831982aSBarry Smith   PetscTruth match;
18273a40ed3dSBarry Smith 
18283a40ed3dSBarry Smith   PetscFunctionBegin;
18294482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18304482741eSBarry Smith   PetscValidCharPointer(type,2);
183182bf6240SBarry Smith 
18326831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
18330f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
183482bf6240SBarry Smith 
183582bf6240SBarry Smith   if (snes->setupcalled) {
1836e1311b90SBarry Smith     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
183782bf6240SBarry Smith     snes->data = 0;
183882bf6240SBarry Smith   }
18397d1a2b2bSBarry Smith 
18409b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
184182bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
184282bf6240SBarry Smith 
1843b9617806SBarry Smith   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
184482bf6240SBarry Smith 
184529bbc08cSBarry Smith   if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type);
18461548b14aSBarry Smith 
1847606d414cSSatish Balay   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
184882bf6240SBarry Smith   snes->data = 0;
18493a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
1850454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
185182bf6240SBarry Smith 
18523a40ed3dSBarry Smith   PetscFunctionReturn(0);
18539b94acceSBarry Smith }
18549b94acceSBarry Smith 
1855a847f771SSatish Balay 
18569b94acceSBarry Smith /* --------------------------------------------------------------------- */
18574a2ae208SSatish Balay #undef __FUNCT__
18584a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
18599b94acceSBarry Smith /*@C
18609b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1861f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
18629b94acceSBarry Smith 
1863fee21e36SBarry Smith    Not Collective
1864fee21e36SBarry Smith 
186536851e7fSLois Curfman McInnes    Level: advanced
186636851e7fSLois Curfman McInnes 
18679b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
18689b94acceSBarry Smith 
18699b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
18709b94acceSBarry Smith @*/
1871cf256101SBarry Smith int SNESRegisterDestroy(void)
18729b94acceSBarry Smith {
187382bf6240SBarry Smith   int ierr;
187482bf6240SBarry Smith 
18753a40ed3dSBarry Smith   PetscFunctionBegin;
187682bf6240SBarry Smith   if (SNESList) {
1877b0a32e0cSBarry Smith     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
187882bf6240SBarry Smith     SNESList = 0;
18799b94acceSBarry Smith   }
18804c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
18813a40ed3dSBarry Smith   PetscFunctionReturn(0);
18829b94acceSBarry Smith }
18839b94acceSBarry Smith 
18844a2ae208SSatish Balay #undef __FUNCT__
18854a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
18869b94acceSBarry Smith /*@C
18879a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
18889b94acceSBarry Smith 
1889c7afd0dbSLois Curfman McInnes    Not Collective
1890c7afd0dbSLois Curfman McInnes 
18919b94acceSBarry Smith    Input Parameter:
18924b0e389bSBarry Smith .  snes - nonlinear solver context
18939b94acceSBarry Smith 
18949b94acceSBarry Smith    Output Parameter:
18953a7fca6bSBarry Smith .  type - SNES method (a character string)
18969b94acceSBarry Smith 
189736851e7fSLois Curfman McInnes    Level: intermediate
189836851e7fSLois Curfman McInnes 
1899454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
19009b94acceSBarry Smith @*/
1901454a90a3SBarry Smith int SNESGetType(SNES snes,SNESType *type)
19029b94acceSBarry Smith {
19033a40ed3dSBarry Smith   PetscFunctionBegin;
19044482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19054482741eSBarry Smith   PetscValidPointer(type,2);
1906454a90a3SBarry Smith   *type = snes->type_name;
19073a40ed3dSBarry Smith   PetscFunctionReturn(0);
19089b94acceSBarry Smith }
19099b94acceSBarry Smith 
19104a2ae208SSatish Balay #undef __FUNCT__
19114a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
19129b94acceSBarry Smith /*@C
19139b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
19149b94acceSBarry Smith    stored.
19159b94acceSBarry Smith 
1916c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1917c7afd0dbSLois Curfman McInnes 
19189b94acceSBarry Smith    Input Parameter:
19199b94acceSBarry Smith .  snes - the SNES context
19209b94acceSBarry Smith 
19219b94acceSBarry Smith    Output Parameter:
19229b94acceSBarry Smith .  x - the solution
19239b94acceSBarry Smith 
192436851e7fSLois Curfman McInnes    Level: advanced
192536851e7fSLois Curfman McInnes 
19269b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
19279b94acceSBarry Smith 
19284b27c08aSLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetSolutionUpdate()
19299b94acceSBarry Smith @*/
19309b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x)
19319b94acceSBarry Smith {
19323a40ed3dSBarry Smith   PetscFunctionBegin;
19334482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19344482741eSBarry Smith   PetscValidPointer(x,2);
19359b94acceSBarry Smith   *x = snes->vec_sol_always;
19363a40ed3dSBarry Smith   PetscFunctionReturn(0);
19379b94acceSBarry Smith }
19389b94acceSBarry Smith 
19394a2ae208SSatish Balay #undef __FUNCT__
19404a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
19419b94acceSBarry Smith /*@C
19429b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
19439b94acceSBarry Smith    stored.
19449b94acceSBarry Smith 
1945c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1946c7afd0dbSLois Curfman McInnes 
19479b94acceSBarry Smith    Input Parameter:
19489b94acceSBarry Smith .  snes - the SNES context
19499b94acceSBarry Smith 
19509b94acceSBarry Smith    Output Parameter:
19519b94acceSBarry Smith .  x - the solution update
19529b94acceSBarry Smith 
195336851e7fSLois Curfman McInnes    Level: advanced
195436851e7fSLois Curfman McInnes 
19559b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
19569b94acceSBarry Smith 
19579b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
19589b94acceSBarry Smith @*/
19599b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x)
19609b94acceSBarry Smith {
19613a40ed3dSBarry Smith   PetscFunctionBegin;
19624482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19634482741eSBarry Smith   PetscValidPointer(x,2);
19649b94acceSBarry Smith   *x = snes->vec_sol_update_always;
19653a40ed3dSBarry Smith   PetscFunctionReturn(0);
19669b94acceSBarry Smith }
19679b94acceSBarry Smith 
19684a2ae208SSatish Balay #undef __FUNCT__
19694a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
19709b94acceSBarry Smith /*@C
19713638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
19729b94acceSBarry Smith 
1973c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1974c7afd0dbSLois Curfman McInnes 
19759b94acceSBarry Smith    Input Parameter:
19769b94acceSBarry Smith .  snes - the SNES context
19779b94acceSBarry Smith 
19789b94acceSBarry Smith    Output Parameter:
19797bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
198000036973SBarry Smith .  ctx - the function context (or PETSC_NULL)
198100036973SBarry Smith -  func - the function (or PETSC_NULL)
19829b94acceSBarry Smith 
198336851e7fSLois Curfman McInnes    Level: advanced
198436851e7fSLois Curfman McInnes 
1985a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
19869b94acceSBarry Smith 
19874b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution()
19889b94acceSBarry Smith @*/
198900036973SBarry Smith int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*))
19909b94acceSBarry Smith {
19913a40ed3dSBarry Smith   PetscFunctionBegin;
19924482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19937bf4e008SBarry Smith   if (r)    *r    = snes->vec_func_always;
19947bf4e008SBarry Smith   if (ctx)  *ctx  = snes->funP;
199500036973SBarry Smith   if (func) *func = snes->computefunction;
19963a40ed3dSBarry Smith   PetscFunctionReturn(0);
19979b94acceSBarry Smith }
19989b94acceSBarry Smith 
19994a2ae208SSatish Balay #undef __FUNCT__
20004a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
20013c7409f5SSatish Balay /*@C
20023c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2003d850072dSLois Curfman McInnes    SNES options in the database.
20043c7409f5SSatish Balay 
2005fee21e36SBarry Smith    Collective on SNES
2006fee21e36SBarry Smith 
2007c7afd0dbSLois Curfman McInnes    Input Parameter:
2008c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2009c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2010c7afd0dbSLois Curfman McInnes 
2011d850072dSLois Curfman McInnes    Notes:
2012a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2013c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2014d850072dSLois Curfman McInnes 
201536851e7fSLois Curfman McInnes    Level: advanced
201636851e7fSLois Curfman McInnes 
20173c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2018a86d99e1SLois Curfman McInnes 
2019a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
20203c7409f5SSatish Balay @*/
20211836bdbcSSatish Balay int SNESSetOptionsPrefix(SNES snes,const char prefix[])
20223c7409f5SSatish Balay {
20233c7409f5SSatish Balay   int ierr;
20243c7409f5SSatish Balay 
20253a40ed3dSBarry Smith   PetscFunctionBegin;
20264482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2027639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
202894b7f48cSBarry Smith   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20293a40ed3dSBarry Smith   PetscFunctionReturn(0);
20303c7409f5SSatish Balay }
20313c7409f5SSatish Balay 
20324a2ae208SSatish Balay #undef __FUNCT__
20334a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
20343c7409f5SSatish Balay /*@C
2035f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2036d850072dSLois Curfman McInnes    SNES options in the database.
20373c7409f5SSatish Balay 
2038fee21e36SBarry Smith    Collective on SNES
2039fee21e36SBarry Smith 
2040c7afd0dbSLois Curfman McInnes    Input Parameters:
2041c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2042c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2043c7afd0dbSLois Curfman McInnes 
2044d850072dSLois Curfman McInnes    Notes:
2045a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2046c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2047d850072dSLois Curfman McInnes 
204836851e7fSLois Curfman McInnes    Level: advanced
204936851e7fSLois Curfman McInnes 
20503c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2051a86d99e1SLois Curfman McInnes 
2052a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20533c7409f5SSatish Balay @*/
20541836bdbcSSatish Balay int SNESAppendOptionsPrefix(SNES snes,const char prefix[])
20553c7409f5SSatish Balay {
20563c7409f5SSatish Balay   int ierr;
20573c7409f5SSatish Balay 
20583a40ed3dSBarry Smith   PetscFunctionBegin;
20594482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2060639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
206194b7f48cSBarry Smith   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20623a40ed3dSBarry Smith   PetscFunctionReturn(0);
20633c7409f5SSatish Balay }
20643c7409f5SSatish Balay 
20654a2ae208SSatish Balay #undef __FUNCT__
20664a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
20679ab63eb5SSatish Balay /*@C
20683c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20693c7409f5SSatish Balay    SNES options in the database.
20703c7409f5SSatish Balay 
2071c7afd0dbSLois Curfman McInnes    Not Collective
2072c7afd0dbSLois Curfman McInnes 
20733c7409f5SSatish Balay    Input Parameter:
20743c7409f5SSatish Balay .  snes - the SNES context
20753c7409f5SSatish Balay 
20763c7409f5SSatish Balay    Output Parameter:
20773c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20783c7409f5SSatish Balay 
20799ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
20809ab63eb5SSatish Balay    sufficient length to hold the prefix.
20819ab63eb5SSatish Balay 
208236851e7fSLois Curfman McInnes    Level: advanced
208336851e7fSLois Curfman McInnes 
20843c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2085a86d99e1SLois Curfman McInnes 
2086a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20873c7409f5SSatish Balay @*/
20881836bdbcSSatish Balay int SNESGetOptionsPrefix(SNES snes,char *prefix[])
20893c7409f5SSatish Balay {
20903c7409f5SSatish Balay   int ierr;
20913c7409f5SSatish Balay 
20923a40ed3dSBarry Smith   PetscFunctionBegin;
20934482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2094639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
20953a40ed3dSBarry Smith   PetscFunctionReturn(0);
20963c7409f5SSatish Balay }
20973c7409f5SSatish Balay 
2098b2002411SLois Curfman McInnes 
20994a2ae208SSatish Balay #undef __FUNCT__
21004a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
21013cea93caSBarry Smith /*@C
21023cea93caSBarry Smith   SNESRegister - See SNESRegisterDynamic()
21033cea93caSBarry Smith 
21047f6c08e0SMatthew Knepley   Level: advanced
21053cea93caSBarry Smith @*/
21061836bdbcSSatish Balay int SNESRegister(const char sname[],const char path[],const char name[],int (*function)(SNES))
2107b2002411SLois Curfman McInnes {
2108b2002411SLois Curfman McInnes   char fullname[256];
2109b2002411SLois Curfman McInnes   int  ierr;
2110b2002411SLois Curfman McInnes 
2111b2002411SLois Curfman McInnes   PetscFunctionBegin;
2112b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2113c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2114b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2115b2002411SLois Curfman McInnes }
2116da9b6338SBarry Smith 
2117da9b6338SBarry Smith #undef __FUNCT__
2118da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin"
2119da9b6338SBarry Smith int SNESTestLocalMin(SNES snes)
2120da9b6338SBarry Smith {
2121da9b6338SBarry Smith   int         ierr,N,i,j;
2122da9b6338SBarry Smith   Vec         u,uh,fh;
2123da9b6338SBarry Smith   PetscScalar value;
2124da9b6338SBarry Smith   PetscReal   norm;
2125da9b6338SBarry Smith 
2126da9b6338SBarry Smith   PetscFunctionBegin;
2127da9b6338SBarry Smith   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2128da9b6338SBarry Smith   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2129da9b6338SBarry Smith   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2130da9b6338SBarry Smith 
2131da9b6338SBarry Smith   /* currently only works for sequential */
2132da9b6338SBarry Smith   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2133da9b6338SBarry Smith   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2134da9b6338SBarry Smith   for (i=0; i<N; i++) {
2135da9b6338SBarry Smith     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2136da9b6338SBarry Smith     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %d\n",i);CHKERRQ(ierr);
2137da9b6338SBarry Smith     for (j=-10; j<11; j++) {
2138ccae9161SBarry Smith       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2139da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
21403ab0aad5SBarry Smith       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2141da9b6338SBarry Smith       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2142da9b6338SBarry Smith       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %d %18.16e\n",j,norm);CHKERRQ(ierr);
2143da9b6338SBarry Smith       value = -value;
2144da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2145da9b6338SBarry Smith     }
2146da9b6338SBarry Smith   }
2147da9b6338SBarry Smith   ierr = VecDestroy(uh);CHKERRQ(ierr);
2148da9b6338SBarry Smith   ierr = VecDestroy(fh);CHKERRQ(ierr);
2149da9b6338SBarry Smith   PetscFunctionReturn(0);
2150da9b6338SBarry Smith }
2151