xref: /petsc/src/snes/interface/snes.c (revision 6849ba73f22fecb8f92ef896a42e4e8bd4cd6965)
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 */
8*6849ba73SBarry Smith PetscCookie SNES_COOKIE = 0;
9*6849ba73SBarry Smith PetscEvent    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 @*/
42dfbe8321SBarry Smith PetscErrorCode SNESView(SNES snes,PetscViewer viewer)
439b94acceSBarry Smith {
449b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
45dfbe8321SBarry Smith   PetscErrorCode ierr;
4694b7f48cSBarry Smith   KSP                 ksp;
47454a90a3SBarry Smith   char                *type;
4832077d6dSBarry 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 
5632077d6dSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr);
57b0a32e0cSBarry Smith   ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr);
5832077d6dSBarry 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;
105*6849ba73SBarry Smith static PetscErrorCode (*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 @*/
121*6849ba73SBarry Smith PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*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 @*/
184dfbe8321SBarry Smith PetscErrorCode SNESSetFromOptions(SNES snes)
1859b94acceSBarry Smith {
18694b7f48cSBarry Smith   KSP                 ksp;
187186905e3SBarry Smith   SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx;
188f1af5d2fSBarry Smith   PetscTruth          flg;
189dfbe8321SBarry Smith   PetscErrorCode ierr;
190dfbe8321SBarry Smith   int  i;
1912fc52814SBarry Smith   const char          *deft;
1922fc52814SBarry Smith   char                type[256];
1939b94acceSBarry Smith 
1943a40ed3dSBarry Smith   PetscFunctionBegin;
1954482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
196ca161407SBarry Smith 
197b0a32e0cSBarry Smith   ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr);
198186905e3SBarry Smith     if (snes->type_name) {
199186905e3SBarry Smith       deft = snes->type_name;
200186905e3SBarry Smith     } else {
2014b27c08aSLois Curfman McInnes       deft = SNESLS;
202d64ed03dSBarry Smith     }
2034bbc92c1SBarry Smith 
204186905e3SBarry Smith     if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
205b0a32e0cSBarry Smith     ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr);
206d64ed03dSBarry Smith     if (flg) {
207186905e3SBarry Smith       ierr = SNESSetType(snes,type);CHKERRQ(ierr);
208186905e3SBarry Smith     } else if (!snes->type_name) {
209186905e3SBarry Smith       ierr = SNESSetType(snes,deft);CHKERRQ(ierr);
210d64ed03dSBarry Smith     }
211909c8a9fSBarry Smith     ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr);
21293c39befSBarry Smith 
21387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr);
21487828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr);
215186905e3SBarry Smith 
21687828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr);
217b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr);
218b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr);
21950ffb88aSMatthew Knepley     ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr);
2205968eb51SBarry Smith     ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr);
2215968eb51SBarry Smith     if (flg) {
2225968eb51SBarry Smith       snes->printreason = PETSC_TRUE;
2235968eb51SBarry Smith     }
224186905e3SBarry Smith 
225b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr);
226186905e3SBarry Smith 
227b0a32e0cSBarry Smith     ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr);
22887828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr);
22987828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr);
23087828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr);
23187828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr);
23287828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr);
23387828ca2SBarry Smith     ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr);
234186905e3SBarry Smith 
235b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr);
23693c39befSBarry Smith     if (flg) {snes->converged = 0;}
237b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr);
2385cd90555SBarry Smith     if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);}
239b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr);
240b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);}
2413a7fca6bSBarry Smith     ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr);
2423a7fca6bSBarry Smith     if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);}
243b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr);
244b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);}
245b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr);
246b8a78c4aSBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);}
247b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr);
2487c922b88SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);}
2495ed2d596SBarry Smith     ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr);
2505ed2d596SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);}
251b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr);
252186905e3SBarry Smith     if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);}
253e24b481bSBarry Smith 
254b0a32e0cSBarry Smith     ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr);
2554b27c08aSLois Curfman McInnes     if (flg) {
256186905e3SBarry Smith       ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr);
257b0a32e0cSBarry Smith       PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n");
2589b94acceSBarry Smith     }
259639f9d9dSBarry Smith 
26076b2cf59SMatthew Knepley     for(i = 0; i < numberofsetfromoptions; i++) {
26176b2cf59SMatthew Knepley       ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr);
26276b2cf59SMatthew Knepley     }
26376b2cf59SMatthew Knepley 
264186905e3SBarry Smith     if (snes->setfromoptions) {
265186905e3SBarry Smith       ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr);
266639f9d9dSBarry Smith     }
267639f9d9dSBarry Smith 
268b0a32e0cSBarry Smith   ierr = PetscOptionsEnd();CHKERRQ(ierr);
2694bbc92c1SBarry Smith 
27094b7f48cSBarry Smith   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
27194b7f48cSBarry Smith   ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr);
27293993e2dSLois Curfman McInnes 
2733a40ed3dSBarry Smith   PetscFunctionReturn(0);
2749b94acceSBarry Smith }
2759b94acceSBarry Smith 
276a847f771SSatish Balay 
2774a2ae208SSatish Balay #undef __FUNCT__
2784a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext"
2799b94acceSBarry Smith /*@
2809b94acceSBarry Smith    SNESSetApplicationContext - Sets the optional user-defined context for
2819b94acceSBarry Smith    the nonlinear solvers.
2829b94acceSBarry Smith 
283fee21e36SBarry Smith    Collective on SNES
284fee21e36SBarry Smith 
285c7afd0dbSLois Curfman McInnes    Input Parameters:
286c7afd0dbSLois Curfman McInnes +  snes - the SNES context
287c7afd0dbSLois Curfman McInnes -  usrP - optional user context
288c7afd0dbSLois Curfman McInnes 
28936851e7fSLois Curfman McInnes    Level: intermediate
29036851e7fSLois Curfman McInnes 
2919b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context
2929b94acceSBarry Smith 
2939b94acceSBarry Smith .seealso: SNESGetApplicationContext()
2949b94acceSBarry Smith @*/
295dfbe8321SBarry Smith PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP)
2969b94acceSBarry Smith {
2973a40ed3dSBarry Smith   PetscFunctionBegin;
2984482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2999b94acceSBarry Smith   snes->user		= usrP;
3003a40ed3dSBarry Smith   PetscFunctionReturn(0);
3019b94acceSBarry Smith }
30274679c65SBarry Smith 
3034a2ae208SSatish Balay #undef __FUNCT__
3044a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext"
3059b94acceSBarry Smith /*@C
3069b94acceSBarry Smith    SNESGetApplicationContext - Gets the user-defined context for the
3079b94acceSBarry Smith    nonlinear solvers.
3089b94acceSBarry Smith 
309c7afd0dbSLois Curfman McInnes    Not Collective
310c7afd0dbSLois Curfman McInnes 
3119b94acceSBarry Smith    Input Parameter:
3129b94acceSBarry Smith .  snes - SNES context
3139b94acceSBarry Smith 
3149b94acceSBarry Smith    Output Parameter:
3159b94acceSBarry Smith .  usrP - user context
3169b94acceSBarry Smith 
31736851e7fSLois Curfman McInnes    Level: intermediate
31836851e7fSLois Curfman McInnes 
3199b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context
3209b94acceSBarry Smith 
3219b94acceSBarry Smith .seealso: SNESSetApplicationContext()
3229b94acceSBarry Smith @*/
323dfbe8321SBarry Smith PetscErrorCode SNESGetApplicationContext(SNES snes,void **usrP)
3249b94acceSBarry Smith {
3253a40ed3dSBarry Smith   PetscFunctionBegin;
3264482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3279b94acceSBarry Smith   *usrP = snes->user;
3283a40ed3dSBarry Smith   PetscFunctionReturn(0);
3299b94acceSBarry Smith }
33074679c65SBarry Smith 
3314a2ae208SSatish Balay #undef __FUNCT__
3324a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber"
3339b94acceSBarry Smith /*@
334c8228a4eSBarry Smith    SNESGetIterationNumber - Gets the number of nonlinear iterations completed
335c8228a4eSBarry Smith    at this time.
3369b94acceSBarry Smith 
337c7afd0dbSLois Curfman McInnes    Not Collective
338c7afd0dbSLois Curfman McInnes 
3399b94acceSBarry Smith    Input Parameter:
3409b94acceSBarry Smith .  snes - SNES context
3419b94acceSBarry Smith 
3429b94acceSBarry Smith    Output Parameter:
3439b94acceSBarry Smith .  iter - iteration number
3449b94acceSBarry Smith 
345c8228a4eSBarry Smith    Notes:
346c8228a4eSBarry Smith    For example, during the computation of iteration 2 this would return 1.
347c8228a4eSBarry Smith 
348c8228a4eSBarry Smith    This is useful for using lagged Jacobians (where one does not recompute the
34908405cd6SLois Curfman McInnes    Jacobian at each SNES iteration). For example, the code
35008405cd6SLois Curfman McInnes .vb
35108405cd6SLois Curfman McInnes       ierr = SNESGetIterationNumber(snes,&it);
35208405cd6SLois Curfman McInnes       if (!(it % 2)) {
35308405cd6SLois Curfman McInnes         [compute Jacobian here]
35408405cd6SLois Curfman McInnes       }
35508405cd6SLois Curfman McInnes .ve
356c8228a4eSBarry Smith    can be used in your ComputeJacobian() function to cause the Jacobian to be
35708405cd6SLois Curfman McInnes    recomputed every second SNES iteration.
358c8228a4eSBarry Smith 
35936851e7fSLois Curfman McInnes    Level: intermediate
36036851e7fSLois Curfman McInnes 
3619b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number
3629b94acceSBarry Smith @*/
363dfbe8321SBarry Smith PetscErrorCode SNESGetIterationNumber(SNES snes,int* iter)
3649b94acceSBarry Smith {
3653a40ed3dSBarry Smith   PetscFunctionBegin;
3664482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3674482741eSBarry Smith   PetscValidIntPointer(iter,2);
3689b94acceSBarry Smith   *iter = snes->iter;
3693a40ed3dSBarry Smith   PetscFunctionReturn(0);
3709b94acceSBarry Smith }
37174679c65SBarry Smith 
3724a2ae208SSatish Balay #undef __FUNCT__
3734a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm"
3749b94acceSBarry Smith /*@
3759b94acceSBarry Smith    SNESGetFunctionNorm - Gets the norm of the current function that was set
3769b94acceSBarry Smith    with SNESSSetFunction().
3779b94acceSBarry Smith 
378c7afd0dbSLois Curfman McInnes    Collective on SNES
379c7afd0dbSLois Curfman McInnes 
3809b94acceSBarry Smith    Input Parameter:
3819b94acceSBarry Smith .  snes - SNES context
3829b94acceSBarry Smith 
3839b94acceSBarry Smith    Output Parameter:
3849b94acceSBarry Smith .  fnorm - 2-norm of function
3859b94acceSBarry Smith 
38636851e7fSLois Curfman McInnes    Level: intermediate
38736851e7fSLois Curfman McInnes 
3889b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm
389a86d99e1SLois Curfman McInnes 
39008405cd6SLois Curfman McInnes .seealso: SNESGetFunction()
3919b94acceSBarry Smith @*/
392dfbe8321SBarry Smith PetscErrorCode SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm)
3939b94acceSBarry Smith {
3943a40ed3dSBarry Smith   PetscFunctionBegin;
3954482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
3964482741eSBarry Smith   PetscValidScalarPointer(fnorm,2);
3979b94acceSBarry Smith   *fnorm = snes->norm;
3983a40ed3dSBarry Smith   PetscFunctionReturn(0);
3999b94acceSBarry Smith }
40074679c65SBarry Smith 
4014a2ae208SSatish Balay #undef __FUNCT__
4024a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps"
4039b94acceSBarry Smith /*@
4049b94acceSBarry Smith    SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps
4059b94acceSBarry Smith    attempted by the nonlinear solver.
4069b94acceSBarry Smith 
407c7afd0dbSLois Curfman McInnes    Not Collective
408c7afd0dbSLois Curfman McInnes 
4099b94acceSBarry Smith    Input Parameter:
4109b94acceSBarry Smith .  snes - SNES context
4119b94acceSBarry Smith 
4129b94acceSBarry Smith    Output Parameter:
4139b94acceSBarry Smith .  nfails - number of unsuccessful steps attempted
4149b94acceSBarry Smith 
415c96a6f78SLois Curfman McInnes    Notes:
416c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
417c96a6f78SLois Curfman McInnes 
41836851e7fSLois Curfman McInnes    Level: intermediate
41936851e7fSLois Curfman McInnes 
4209b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps
4219b94acceSBarry Smith @*/
422dfbe8321SBarry Smith PetscErrorCode SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails)
4239b94acceSBarry Smith {
4243a40ed3dSBarry Smith   PetscFunctionBegin;
4254482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4264482741eSBarry Smith   PetscValidIntPointer(nfails,2);
42750ffb88aSMatthew Knepley   *nfails = snes->numFailures;
42850ffb88aSMatthew Knepley   PetscFunctionReturn(0);
42950ffb88aSMatthew Knepley }
43050ffb88aSMatthew Knepley 
43150ffb88aSMatthew Knepley #undef __FUNCT__
43250ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps"
43350ffb88aSMatthew Knepley /*@
43450ffb88aSMatthew Knepley    SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps
43550ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
43650ffb88aSMatthew Knepley 
43750ffb88aSMatthew Knepley    Not Collective
43850ffb88aSMatthew Knepley 
43950ffb88aSMatthew Knepley    Input Parameters:
44050ffb88aSMatthew Knepley +  snes     - SNES context
44150ffb88aSMatthew Knepley -  maxFails - maximum of unsuccessful steps
44250ffb88aSMatthew Knepley 
44350ffb88aSMatthew Knepley    Level: intermediate
44450ffb88aSMatthew Knepley 
44550ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps
44650ffb88aSMatthew Knepley @*/
447dfbe8321SBarry Smith PetscErrorCode SNESSetMaximumUnsuccessfulSteps(SNES snes, int maxFails)
44850ffb88aSMatthew Knepley {
44950ffb88aSMatthew Knepley   PetscFunctionBegin;
4504482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
45150ffb88aSMatthew Knepley   snes->maxFailures = maxFails;
45250ffb88aSMatthew Knepley   PetscFunctionReturn(0);
45350ffb88aSMatthew Knepley }
45450ffb88aSMatthew Knepley 
45550ffb88aSMatthew Knepley #undef __FUNCT__
45650ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps"
45750ffb88aSMatthew Knepley /*@
45850ffb88aSMatthew Knepley    SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps
45950ffb88aSMatthew Knepley    attempted by the nonlinear solver before it gives up.
46050ffb88aSMatthew Knepley 
46150ffb88aSMatthew Knepley    Not Collective
46250ffb88aSMatthew Knepley 
46350ffb88aSMatthew Knepley    Input Parameter:
46450ffb88aSMatthew Knepley .  snes     - SNES context
46550ffb88aSMatthew Knepley 
46650ffb88aSMatthew Knepley    Output Parameter:
46750ffb88aSMatthew Knepley .  maxFails - maximum of unsuccessful steps
46850ffb88aSMatthew Knepley 
46950ffb88aSMatthew Knepley    Level: intermediate
47050ffb88aSMatthew Knepley 
47150ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps
47250ffb88aSMatthew Knepley @*/
473dfbe8321SBarry Smith PetscErrorCode SNESGetMaximumUnsuccessfulSteps(SNES snes, int *maxFails)
47450ffb88aSMatthew Knepley {
47550ffb88aSMatthew Knepley   PetscFunctionBegin;
4764482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
4774482741eSBarry Smith   PetscValidIntPointer(maxFails,2);
47850ffb88aSMatthew Knepley   *maxFails = snes->maxFailures;
4793a40ed3dSBarry Smith   PetscFunctionReturn(0);
4809b94acceSBarry Smith }
481a847f771SSatish Balay 
4824a2ae208SSatish Balay #undef __FUNCT__
4834a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations"
484c96a6f78SLois Curfman McInnes /*@
485c96a6f78SLois Curfman McInnes    SNESGetNumberLinearIterations - Gets the total number of linear iterations
486c96a6f78SLois Curfman McInnes    used by the nonlinear solver.
487c96a6f78SLois Curfman McInnes 
488c7afd0dbSLois Curfman McInnes    Not Collective
489c7afd0dbSLois Curfman McInnes 
490c96a6f78SLois Curfman McInnes    Input Parameter:
491c96a6f78SLois Curfman McInnes .  snes - SNES context
492c96a6f78SLois Curfman McInnes 
493c96a6f78SLois Curfman McInnes    Output Parameter:
494c96a6f78SLois Curfman McInnes .  lits - number of linear iterations
495c96a6f78SLois Curfman McInnes 
496c96a6f78SLois Curfman McInnes    Notes:
497c96a6f78SLois Curfman McInnes    This counter is reset to zero for each successive call to SNESSolve().
498c96a6f78SLois Curfman McInnes 
49936851e7fSLois Curfman McInnes    Level: intermediate
50036851e7fSLois Curfman McInnes 
501c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations
502c96a6f78SLois Curfman McInnes @*/
503dfbe8321SBarry Smith PetscErrorCode SNESGetNumberLinearIterations(SNES snes,int* lits)
504c96a6f78SLois Curfman McInnes {
5053a40ed3dSBarry Smith   PetscFunctionBegin;
5064482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5074482741eSBarry Smith   PetscValidIntPointer(lits,2);
508c96a6f78SLois Curfman McInnes   *lits = snes->linear_its;
5093a40ed3dSBarry Smith   PetscFunctionReturn(0);
510c96a6f78SLois Curfman McInnes }
511c96a6f78SLois Curfman McInnes 
5124a2ae208SSatish Balay #undef __FUNCT__
51394b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP"
5149b94acceSBarry Smith /*@C
51594b7f48cSBarry Smith    SNESGetKSP - Returns the KSP context for a SNES solver.
5169b94acceSBarry Smith 
51794b7f48cSBarry Smith    Not Collective, but if SNES object is parallel, then KSP object is parallel
518c7afd0dbSLois Curfman McInnes 
5199b94acceSBarry Smith    Input Parameter:
5209b94acceSBarry Smith .  snes - the SNES context
5219b94acceSBarry Smith 
5229b94acceSBarry Smith    Output Parameter:
52394b7f48cSBarry Smith .  ksp - the KSP context
5249b94acceSBarry Smith 
5259b94acceSBarry Smith    Notes:
52694b7f48cSBarry Smith    The user can then directly manipulate the KSP context to set various
5279b94acceSBarry Smith    options, etc.  Likewise, the user can then extract and manipulate the
5289b94acceSBarry Smith    KSP and PC contexts as well.
5299b94acceSBarry Smith 
53036851e7fSLois Curfman McInnes    Level: beginner
53136851e7fSLois Curfman McInnes 
53294b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context
5339b94acceSBarry Smith 
53494b7f48cSBarry Smith .seealso: KSPGetPC()
5359b94acceSBarry Smith @*/
536dfbe8321SBarry Smith PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp)
5379b94acceSBarry Smith {
5383a40ed3dSBarry Smith   PetscFunctionBegin;
5394482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
5404482741eSBarry Smith   PetscValidPointer(ksp,2);
54194b7f48cSBarry Smith   *ksp = snes->ksp;
5423a40ed3dSBarry Smith   PetscFunctionReturn(0);
5439b94acceSBarry Smith }
54482bf6240SBarry Smith 
5454a2ae208SSatish Balay #undef __FUNCT__
5464a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc"
547*6849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj)
548e24b481bSBarry Smith {
549aa482453SBarry Smith #if defined(PETSC_HAVE_AMS)
550454a90a3SBarry Smith   SNES          v = (SNES) obj;
551dfbe8321SBarry Smith   PetscErrorCode ierr;
55243d6d2cbSBarry Smith #endif
553e24b481bSBarry Smith 
554e24b481bSBarry Smith   PetscFunctionBegin;
555e24b481bSBarry Smith 
55643d6d2cbSBarry Smith #if defined(PETSC_HAVE_AMS)
557e24b481bSBarry Smith   /* if it is already published then return */
558e24b481bSBarry Smith   if (v->amem >=0) PetscFunctionReturn(0);
559e24b481bSBarry Smith 
560454a90a3SBarry Smith   ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr);
561e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ,
562e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
563e24b481bSBarry Smith   ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ,
564e24b481bSBarry Smith                                 AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr);
565454a90a3SBarry Smith   ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr);
566433994e6SBarry Smith #endif
567e24b481bSBarry Smith   PetscFunctionReturn(0);
568e24b481bSBarry Smith }
569e24b481bSBarry Smith 
5709b94acceSBarry Smith /* -----------------------------------------------------------*/
5714a2ae208SSatish Balay #undef __FUNCT__
5724a2ae208SSatish Balay #define __FUNCT__ "SNESCreate"
5739b94acceSBarry Smith /*@C
5749b94acceSBarry Smith    SNESCreate - Creates a nonlinear solver context.
5759b94acceSBarry Smith 
576c7afd0dbSLois Curfman McInnes    Collective on MPI_Comm
577c7afd0dbSLois Curfman McInnes 
578c7afd0dbSLois Curfman McInnes    Input Parameters:
579c7afd0dbSLois Curfman McInnes +  comm - MPI communicator
5809b94acceSBarry Smith 
5819b94acceSBarry Smith    Output Parameter:
5829b94acceSBarry Smith .  outsnes - the new SNES context
5839b94acceSBarry Smith 
584c7afd0dbSLois Curfman McInnes    Options Database Keys:
585c7afd0dbSLois Curfman McInnes +   -snes_mf - Activates default matrix-free Jacobian-vector products,
586c7afd0dbSLois Curfman McInnes                and no preconditioning matrix
587c7afd0dbSLois Curfman McInnes .   -snes_mf_operator - Activates default matrix-free Jacobian-vector
588c7afd0dbSLois Curfman McInnes                products, and a user-provided preconditioning matrix
589c7afd0dbSLois Curfman McInnes                as set by SNESSetJacobian()
590c7afd0dbSLois Curfman McInnes -   -snes_fd - Uses (slow!) finite differences to compute Jacobian
591c1f60f51SBarry Smith 
59236851e7fSLois Curfman McInnes    Level: beginner
59336851e7fSLois Curfman McInnes 
5949b94acceSBarry Smith .keywords: SNES, nonlinear, create, context
5959b94acceSBarry Smith 
5964b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES
5979b94acceSBarry Smith @*/
598dfbe8321SBarry Smith PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes)
5999b94acceSBarry Smith {
600dfbe8321SBarry Smith   PetscErrorCode ierr;
6019b94acceSBarry Smith   SNES                snes;
6029b94acceSBarry Smith   SNES_KSP_EW_ConvCtx *kctx;
60337fcc0dbSBarry Smith 
6043a40ed3dSBarry Smith   PetscFunctionBegin;
6054482741eSBarry Smith   PetscValidPointer(outsnes,1);
6068ba1e511SMatthew Knepley   *outsnes = PETSC_NULL;
6078ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES
6088ba1e511SMatthew Knepley   ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr);
6098ba1e511SMatthew Knepley #endif
6108ba1e511SMatthew Knepley 
6113f1db9ecSBarry Smith   PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);
612b0a32e0cSBarry Smith   PetscLogObjectCreate(snes);
613e24b481bSBarry Smith   snes->bops->publish     = SNESPublish_Petsc;
6149b94acceSBarry Smith   snes->max_its           = 50;
6159750a799SBarry Smith   snes->max_funcs	  = 10000;
6169b94acceSBarry Smith   snes->norm		  = 0.0;
617b4874afaSBarry Smith   snes->rtol		  = 1.e-8;
618b4874afaSBarry Smith   snes->ttol              = 0.0;
619b4874afaSBarry Smith   snes->atol		  = 1.e-50;
6209b94acceSBarry Smith   snes->xtol		  = 1.e-8;
6214b27c08aSLois Curfman McInnes   snes->deltatol	  = 1.e-12;
6229b94acceSBarry Smith   snes->nfuncs            = 0;
62350ffb88aSMatthew Knepley   snes->numFailures       = 0;
62450ffb88aSMatthew Knepley   snes->maxFailures       = 1;
6257a00f4a9SLois Curfman McInnes   snes->linear_its        = 0;
626639f9d9dSBarry Smith   snes->numbermonitors    = 0;
6279b94acceSBarry Smith   snes->data              = 0;
6289b94acceSBarry Smith   snes->view              = 0;
62982bf6240SBarry Smith   snes->setupcalled       = 0;
630186905e3SBarry Smith   snes->ksp_ewconv        = PETSC_FALSE;
6316f24a144SLois Curfman McInnes   snes->vwork             = 0;
6326f24a144SLois Curfman McInnes   snes->nwork             = 0;
633758f92a0SBarry Smith   snes->conv_hist_len     = 0;
634758f92a0SBarry Smith   snes->conv_hist_max     = 0;
635758f92a0SBarry Smith   snes->conv_hist         = PETSC_NULL;
636758f92a0SBarry Smith   snes->conv_hist_its     = PETSC_NULL;
637758f92a0SBarry Smith   snes->conv_hist_reset   = PETSC_TRUE;
638184914b5SBarry Smith   snes->reason            = SNES_CONVERGED_ITERATING;
6399b94acceSBarry Smith 
6409b94acceSBarry Smith   /* Create context to compute Eisenstat-Walker relative tolerance for KSP */
641b0a32e0cSBarry Smith   ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr);
642b0a32e0cSBarry Smith   PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));
6439b94acceSBarry Smith   snes->kspconvctx  = (void*)kctx;
6449b94acceSBarry Smith   kctx->version     = 2;
6459b94acceSBarry Smith   kctx->rtol_0      = .3; /* Eisenstat and Walker suggest rtol_0=.5, but
6469b94acceSBarry Smith                              this was too large for some test cases */
6479b94acceSBarry Smith   kctx->rtol_last   = 0;
6489b94acceSBarry Smith   kctx->rtol_max    = .9;
6499b94acceSBarry Smith   kctx->gamma       = 1.0;
6509b94acceSBarry Smith   kctx->alpha2      = .5*(1.0 + sqrt(5.0));
6519b94acceSBarry Smith   kctx->alpha       = kctx->alpha2;
6529b94acceSBarry Smith   kctx->threshold   = .1;
6539b94acceSBarry Smith   kctx->lresid_last = 0;
6549b94acceSBarry Smith   kctx->norm_last   = 0;
6559b94acceSBarry Smith 
65694b7f48cSBarry Smith   ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr);
65794b7f48cSBarry Smith   PetscLogObjectParent(snes,snes->ksp)
6585334005bSBarry Smith 
6599b94acceSBarry Smith   *outsnes = snes;
66000036973SBarry Smith   ierr = PetscPublishAll(snes);CHKERRQ(ierr);
6613a40ed3dSBarry Smith   PetscFunctionReturn(0);
6629b94acceSBarry Smith }
6639b94acceSBarry Smith 
6644a2ae208SSatish Balay #undef __FUNCT__
6654a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction"
6669b94acceSBarry Smith /*@C
6679b94acceSBarry Smith    SNESSetFunction - Sets the function evaluation routine and function
6689b94acceSBarry Smith    vector for use by the SNES routines in solving systems of nonlinear
6699b94acceSBarry Smith    equations.
6709b94acceSBarry Smith 
671fee21e36SBarry Smith    Collective on SNES
672fee21e36SBarry Smith 
673c7afd0dbSLois Curfman McInnes    Input Parameters:
674c7afd0dbSLois Curfman McInnes +  snes - the SNES context
675c7afd0dbSLois Curfman McInnes .  func - function evaluation routine
676c7afd0dbSLois Curfman McInnes .  r - vector to store function value
677c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
678c7afd0dbSLois Curfman McInnes          function evaluation routine (may be PETSC_NULL)
6799b94acceSBarry Smith 
680c7afd0dbSLois Curfman McInnes    Calling sequence of func:
6818d76a1e5SLois Curfman McInnes $    func (SNES snes,Vec x,Vec f,void *ctx);
682c7afd0dbSLois Curfman McInnes 
683313e4042SLois Curfman McInnes .  f - function vector
684c7afd0dbSLois Curfman McInnes -  ctx - optional user-defined function context
6859b94acceSBarry Smith 
6869b94acceSBarry Smith    Notes:
6879b94acceSBarry Smith    The Newton-like methods typically solve linear systems of the form
6889b94acceSBarry Smith $      f'(x) x = -f(x),
689c7afd0dbSLois Curfman McInnes    where f'(x) denotes the Jacobian matrix and f(x) is the function.
6909b94acceSBarry Smith 
69136851e7fSLois Curfman McInnes    Level: beginner
69236851e7fSLois Curfman McInnes 
6939b94acceSBarry Smith .keywords: SNES, nonlinear, set, function
6949b94acceSBarry Smith 
695a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian()
6969b94acceSBarry Smith @*/
697*6849ba73SBarry Smith PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx)
6989b94acceSBarry Smith {
6993a40ed3dSBarry Smith   PetscFunctionBegin;
7004482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7014482741eSBarry Smith   PetscValidHeaderSpecific(r,VEC_COOKIE,2);
702c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,r,2);
703184914b5SBarry Smith 
7049b94acceSBarry Smith   snes->computefunction     = func;
7059b94acceSBarry Smith   snes->vec_func            = snes->vec_func_always = r;
7069b94acceSBarry Smith   snes->funP                = ctx;
7073a40ed3dSBarry Smith   PetscFunctionReturn(0);
7089b94acceSBarry Smith }
7099b94acceSBarry Smith 
7103ab0aad5SBarry Smith /* --------------------------------------------------------------- */
7113ab0aad5SBarry Smith #undef __FUNCT__
7123ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs"
7133ab0aad5SBarry Smith /*@C
7143ab0aad5SBarry Smith    SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set
7153ab0aad5SBarry Smith    it assumes a zero right hand side.
7163ab0aad5SBarry Smith 
7173ab0aad5SBarry Smith    Collective on SNES
7183ab0aad5SBarry Smith 
7193ab0aad5SBarry Smith    Input Parameters:
7203ab0aad5SBarry Smith +  snes - the SNES context
7213ab0aad5SBarry Smith -  rhs - the right hand side vector or PETSC_NULL for a zero right hand side
7223ab0aad5SBarry Smith 
7233ab0aad5SBarry Smith    Level: intermediate
7243ab0aad5SBarry Smith 
7253ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side
7263ab0aad5SBarry Smith 
7273ab0aad5SBarry Smith .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction()
7283ab0aad5SBarry Smith @*/
729dfbe8321SBarry Smith PetscErrorCode SNESSetRhs(SNES snes,Vec rhs)
7303ab0aad5SBarry Smith {
731dfbe8321SBarry Smith   PetscErrorCode ierr;
7323ab0aad5SBarry Smith 
7333ab0aad5SBarry Smith   PetscFunctionBegin;
7343ab0aad5SBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7353ab0aad5SBarry Smith   if (rhs) {
7363ab0aad5SBarry Smith     PetscValidHeaderSpecific(rhs,VEC_COOKIE,2);
7373ab0aad5SBarry Smith     PetscCheckSameComm(snes,1,rhs,2);
7383ab0aad5SBarry Smith     ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr);
7393ab0aad5SBarry Smith   }
7403ab0aad5SBarry Smith   if (snes->afine) {
7413ab0aad5SBarry Smith     ierr = VecDestroy(snes->afine);CHKERRQ(ierr);
7423ab0aad5SBarry Smith   }
7433ab0aad5SBarry Smith   snes->afine = rhs;
7443ab0aad5SBarry Smith   PetscFunctionReturn(0);
7453ab0aad5SBarry Smith }
7463ab0aad5SBarry Smith 
7474a2ae208SSatish Balay #undef __FUNCT__
7484a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction"
7499b94acceSBarry Smith /*@
75036851e7fSLois Curfman McInnes    SNESComputeFunction - Calls the function that has been set with
7519b94acceSBarry Smith                          SNESSetFunction().
7529b94acceSBarry Smith 
753c7afd0dbSLois Curfman McInnes    Collective on SNES
754c7afd0dbSLois Curfman McInnes 
7559b94acceSBarry Smith    Input Parameters:
756c7afd0dbSLois Curfman McInnes +  snes - the SNES context
757c7afd0dbSLois Curfman McInnes -  x - input vector
7589b94acceSBarry Smith 
7599b94acceSBarry Smith    Output Parameter:
7603638b69dSLois Curfman McInnes .  y - function vector, as set by SNESSetFunction()
7619b94acceSBarry Smith 
7621bffabb2SLois Curfman McInnes    Notes:
76336851e7fSLois Curfman McInnes    SNESComputeFunction() is typically used within nonlinear solvers
76436851e7fSLois Curfman McInnes    implementations, so most users would not generally call this routine
76536851e7fSLois Curfman McInnes    themselves.
76636851e7fSLois Curfman McInnes 
76736851e7fSLois Curfman McInnes    Level: developer
76836851e7fSLois Curfman McInnes 
7699b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function
7709b94acceSBarry Smith 
771a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction()
7729b94acceSBarry Smith @*/
773dfbe8321SBarry Smith PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y)
7749b94acceSBarry Smith {
775dfbe8321SBarry Smith   PetscErrorCode ierr;
7769b94acceSBarry Smith 
7773a40ed3dSBarry Smith   PetscFunctionBegin;
7784482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
7794482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
7804482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,3);
781c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
782c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,3);
783184914b5SBarry Smith 
784d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
785d64ed03dSBarry Smith   PetscStackPush("SNES user function");
7869b94acceSBarry Smith   ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr);
787d64ed03dSBarry Smith   PetscStackPop;
7883ab0aad5SBarry Smith   if (snes->afine) {
7893ab0aad5SBarry Smith     PetscScalar mone = -1.0;
7903ab0aad5SBarry Smith     ierr = VecAXPY(&mone,snes->afine,y);CHKERRQ(ierr);
7913ab0aad5SBarry Smith   }
792ae3c334cSLois Curfman McInnes   snes->nfuncs++;
793d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr);
7943a40ed3dSBarry Smith   PetscFunctionReturn(0);
7959b94acceSBarry Smith }
7969b94acceSBarry Smith 
7974a2ae208SSatish Balay #undef __FUNCT__
7984a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian"
79962fef451SLois Curfman McInnes /*@
80062fef451SLois Curfman McInnes    SNESComputeJacobian - Computes the Jacobian matrix that has been
80162fef451SLois Curfman McInnes    set with SNESSetJacobian().
80262fef451SLois Curfman McInnes 
803c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
804c7afd0dbSLois Curfman McInnes 
80562fef451SLois Curfman McInnes    Input Parameters:
806c7afd0dbSLois Curfman McInnes +  snes - the SNES context
807c7afd0dbSLois Curfman McInnes -  x - input vector
80862fef451SLois Curfman McInnes 
80962fef451SLois Curfman McInnes    Output Parameters:
810c7afd0dbSLois Curfman McInnes +  A - Jacobian matrix
81162fef451SLois Curfman McInnes .  B - optional preconditioning matrix
812c7afd0dbSLois Curfman McInnes -  flag - flag indicating matrix structure
813fee21e36SBarry Smith 
81462fef451SLois Curfman McInnes    Notes:
81562fef451SLois Curfman McInnes    Most users should not need to explicitly call this routine, as it
81662fef451SLois Curfman McInnes    is used internally within the nonlinear solvers.
81762fef451SLois Curfman McInnes 
81894b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the
819dc5a77f8SLois Curfman McInnes    flag parameter.
82062fef451SLois Curfman McInnes 
82136851e7fSLois Curfman McInnes    Level: developer
82236851e7fSLois Curfman McInnes 
82362fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix
82462fef451SLois Curfman McInnes 
82594b7f48cSBarry Smith .seealso:  SNESSetJacobian(), KSPSetOperators()
82662fef451SLois Curfman McInnes @*/
827dfbe8321SBarry Smith PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg)
8289b94acceSBarry Smith {
829dfbe8321SBarry Smith   PetscErrorCode ierr;
8303a40ed3dSBarry Smith 
8313a40ed3dSBarry Smith   PetscFunctionBegin;
8324482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8334482741eSBarry Smith   PetscValidHeaderSpecific(X,VEC_COOKIE,2);
8344482741eSBarry Smith   PetscValidPointer(flg,5);
835c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,X,2);
8363a40ed3dSBarry Smith   if (!snes->computejacobian) PetscFunctionReturn(0);
837d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
838c4fc05e7SBarry Smith   *flg = DIFFERENT_NONZERO_PATTERN;
839d64ed03dSBarry Smith   PetscStackPush("SNES user Jacobian function");
8409b94acceSBarry Smith   ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr);
841d64ed03dSBarry Smith   PetscStackPop;
842d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr);
8436d84be18SBarry Smith   /* make sure user returned a correct Jacobian and preconditioner */
8444482741eSBarry Smith   PetscValidHeaderSpecific(*A,MAT_COOKIE,3);
8454482741eSBarry Smith   PetscValidHeaderSpecific(*B,MAT_COOKIE,4);
8463a40ed3dSBarry Smith   PetscFunctionReturn(0);
8479b94acceSBarry Smith }
8489b94acceSBarry Smith 
8494a2ae208SSatish Balay #undef __FUNCT__
8504a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian"
8519b94acceSBarry Smith /*@C
8529b94acceSBarry Smith    SNESSetJacobian - Sets the function to compute Jacobian as well as the
853044dda88SLois Curfman McInnes    location to store the matrix.
8549b94acceSBarry Smith 
855c7afd0dbSLois Curfman McInnes    Collective on SNES and Mat
856c7afd0dbSLois Curfman McInnes 
8579b94acceSBarry Smith    Input Parameters:
858c7afd0dbSLois Curfman McInnes +  snes - the SNES context
8599b94acceSBarry Smith .  A - Jacobian matrix
8609b94acceSBarry Smith .  B - preconditioner matrix (usually same as the Jacobian)
8619b94acceSBarry Smith .  func - Jacobian evaluation routine
862c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined context for private data for the
8632cd2dfdcSLois Curfman McInnes          Jacobian evaluation routine (may be PETSC_NULL)
8649b94acceSBarry Smith 
8659b94acceSBarry Smith    Calling sequence of func:
8668d76a1e5SLois Curfman McInnes $     func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx);
8679b94acceSBarry Smith 
868c7afd0dbSLois Curfman McInnes +  x - input vector
8699b94acceSBarry Smith .  A - Jacobian matrix
8709b94acceSBarry Smith .  B - preconditioner matrix, usually the same as A
871ac21db08SLois Curfman McInnes .  flag - flag indicating information about the preconditioner matrix
87294b7f48cSBarry Smith    structure (same as flag in KSPSetOperators())
873c7afd0dbSLois Curfman McInnes -  ctx - [optional] user-defined Jacobian context
8749b94acceSBarry Smith 
8759b94acceSBarry Smith    Notes:
87694b7f48cSBarry Smith    See KSPSetOperators() for important information about setting the flag
8772cd2dfdcSLois Curfman McInnes    output parameter in the routine func().  Be sure to read this information!
878ac21db08SLois Curfman McInnes 
879ac21db08SLois Curfman McInnes    The routine func() takes Mat * as the matrix arguments rather than Mat.
8809b94acceSBarry Smith    This allows the Jacobian evaluation routine to replace A and/or B with a
8819b94acceSBarry Smith    completely new new matrix structure (not just different matrix elements)
8829b94acceSBarry Smith    when appropriate, for instance, if the nonzero structure is changing
8839b94acceSBarry Smith    throughout the global iterations.
8849b94acceSBarry Smith 
88536851e7fSLois Curfman McInnes    Level: beginner
88636851e7fSLois Curfman McInnes 
8879b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix
8889b94acceSBarry Smith 
88994b7f48cSBarry Smith .seealso: KSPSetOperators(), SNESSetFunction()
8909b94acceSBarry Smith @*/
891*6849ba73SBarry Smith PetscErrorCode SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx)
8929b94acceSBarry Smith {
893dfbe8321SBarry Smith   PetscErrorCode ierr;
8943a7fca6bSBarry Smith 
8953a40ed3dSBarry Smith   PetscFunctionBegin;
8964482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
8974482741eSBarry Smith   if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2);
8984482741eSBarry Smith   if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3);
899c9780b6fSBarry Smith   if (A) PetscCheckSameComm(snes,1,A,2);
900c9780b6fSBarry Smith   if (B) PetscCheckSameComm(snes,1,B,2);
9013a7fca6bSBarry Smith   if (func) snes->computejacobian = func;
9023a7fca6bSBarry Smith   if (ctx)  snes->jacP            = ctx;
9033a7fca6bSBarry Smith   if (A) {
9043a7fca6bSBarry Smith     if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
9059b94acceSBarry Smith     snes->jacobian = A;
9063a7fca6bSBarry Smith     ierr           = PetscObjectReference((PetscObject)A);CHKERRQ(ierr);
9073a7fca6bSBarry Smith   }
9083a7fca6bSBarry Smith   if (B) {
9093a7fca6bSBarry Smith     if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
9109b94acceSBarry Smith     snes->jacobian_pre = B;
9113a7fca6bSBarry Smith     ierr               = PetscObjectReference((PetscObject)B);CHKERRQ(ierr);
9123a7fca6bSBarry Smith   }
9133a40ed3dSBarry Smith   PetscFunctionReturn(0);
9149b94acceSBarry Smith }
91562fef451SLois Curfman McInnes 
9164a2ae208SSatish Balay #undef __FUNCT__
9174a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian"
918c2aafc4cSSatish Balay /*@C
919b4fd4287SBarry Smith    SNESGetJacobian - Returns the Jacobian matrix and optionally the user
920b4fd4287SBarry Smith    provided context for evaluating the Jacobian.
921b4fd4287SBarry Smith 
922c7afd0dbSLois Curfman McInnes    Not Collective, but Mat object will be parallel if SNES object is
923c7afd0dbSLois Curfman McInnes 
924b4fd4287SBarry Smith    Input Parameter:
925b4fd4287SBarry Smith .  snes - the nonlinear solver context
926b4fd4287SBarry Smith 
927b4fd4287SBarry Smith    Output Parameters:
928c7afd0dbSLois Curfman McInnes +  A - location to stash Jacobian matrix (or PETSC_NULL)
929b4fd4287SBarry Smith .  B - location to stash preconditioner matrix (or PETSC_NULL)
93000036973SBarry Smith .  ctx - location to stash Jacobian ctx (or PETSC_NULL)
93100036973SBarry Smith -  func - location to put Jacobian function (or PETSC_NULL)
932fee21e36SBarry Smith 
93336851e7fSLois Curfman McInnes    Level: advanced
93436851e7fSLois Curfman McInnes 
935b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian()
936b4fd4287SBarry Smith @*/
937*6849ba73SBarry Smith PetscErrorCode SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*))
938b4fd4287SBarry Smith {
9393a40ed3dSBarry Smith   PetscFunctionBegin;
9404482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
941b4fd4287SBarry Smith   if (A)    *A    = snes->jacobian;
942b4fd4287SBarry Smith   if (B)    *B    = snes->jacobian_pre;
943b4fd4287SBarry Smith   if (ctx)  *ctx  = snes->jacP;
94400036973SBarry Smith   if (func) *func = snes->computejacobian;
9453a40ed3dSBarry Smith   PetscFunctionReturn(0);
946b4fd4287SBarry Smith }
947b4fd4287SBarry Smith 
9489b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */
949dfbe8321SBarry Smith EXTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*);
9509b94acceSBarry Smith 
9514a2ae208SSatish Balay #undef __FUNCT__
9524a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp"
9539b94acceSBarry Smith /*@
9549b94acceSBarry Smith    SNESSetUp - Sets up the internal data structures for the later use
955272ac6f2SLois Curfman McInnes    of a nonlinear solver.
9569b94acceSBarry Smith 
957fee21e36SBarry Smith    Collective on SNES
958fee21e36SBarry Smith 
959c7afd0dbSLois Curfman McInnes    Input Parameters:
960c7afd0dbSLois Curfman McInnes +  snes - the SNES context
961c7afd0dbSLois Curfman McInnes -  x - the solution vector
962c7afd0dbSLois Curfman McInnes 
963272ac6f2SLois Curfman McInnes    Notes:
964272ac6f2SLois Curfman McInnes    For basic use of the SNES solvers the user need not explicitly call
965272ac6f2SLois Curfman McInnes    SNESSetUp(), since these actions will automatically occur during
966272ac6f2SLois Curfman McInnes    the call to SNESSolve().  However, if one wishes to control this
967272ac6f2SLois Curfman McInnes    phase separately, SNESSetUp() should be called after SNESCreate()
968272ac6f2SLois Curfman McInnes    and optional routines of the form SNESSetXXX(), but before SNESSolve().
969272ac6f2SLois Curfman McInnes 
97036851e7fSLois Curfman McInnes    Level: advanced
97136851e7fSLois Curfman McInnes 
9729b94acceSBarry Smith .keywords: SNES, nonlinear, setup
9739b94acceSBarry Smith 
9749b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy()
9759b94acceSBarry Smith @*/
976dfbe8321SBarry Smith PetscErrorCode SNESSetUp(SNES snes,Vec x)
9779b94acceSBarry Smith {
978dfbe8321SBarry Smith   PetscErrorCode ierr;
9794b27c08aSLois Curfman McInnes   PetscTruth flg, iseqtr;
9803a40ed3dSBarry Smith 
9813a40ed3dSBarry Smith   PetscFunctionBegin;
9824482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
9834482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
984c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
9858ddd3da0SLois Curfman McInnes   snes->vec_sol = snes->vec_sol_always = x;
9869b94acceSBarry Smith 
987b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr);
988c1f60f51SBarry Smith   /*
989c1f60f51SBarry Smith       This version replaces the user provided Jacobian matrix with a
990dfa02198SLois Curfman McInnes       matrix-free version but still employs the user-provided preconditioner matrix
991c1f60f51SBarry Smith   */
992c1f60f51SBarry Smith   if (flg) {
993c1f60f51SBarry Smith     Mat J;
9945a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
9955a655dc6SBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
9963a7fca6bSBarry Smith     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n");
9973a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
9983a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
999c1f60f51SBarry Smith   }
100045fc7adcSBarry Smith 
1001b4014bb2SMatthew Knepley #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE)
100245fc7adcSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr);
100345fc7adcSBarry Smith   if (flg) {
100445fc7adcSBarry Smith     Mat J;
100545fc7adcSBarry Smith     ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr);
100645fc7adcSBarry Smith     ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr);
100745fc7adcSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
100845fc7adcSBarry Smith   }
100932a4b47aSMatthew Knepley #endif
101045fc7adcSBarry Smith 
1011b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr);
1012c1f60f51SBarry Smith   /*
1013dfa02198SLois Curfman McInnes       This version replaces both the user-provided Jacobian and the user-
1014c1f60f51SBarry Smith       provided preconditioner matrix with the default matrix free version.
1015c1f60f51SBarry Smith    */
101631615d04SLois Curfman McInnes   if (flg) {
1017272ac6f2SLois Curfman McInnes     Mat  J;
1018b5d62d44SBarry Smith     KSP ksp;
101994b7f48cSBarry Smith     PC   pc;
1020f3ef73ceSBarry Smith 
10215a655dc6SBarry Smith     ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr);
10223a7fca6bSBarry Smith     ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr);
102393b98531SBarry Smith     PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n");
10243a7fca6bSBarry Smith     ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr);
10253a7fca6bSBarry Smith     ierr = MatDestroy(J);CHKERRQ(ierr);
10263a7fca6bSBarry Smith 
1027f3ef73ceSBarry Smith     /* force no preconditioner */
102894b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1029b5d62d44SBarry Smith     ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
1030a9815358SBarry Smith     ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr);
1031a9815358SBarry Smith     if (!flg) {
1032f3ef73ceSBarry Smith       ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr);
1033272ac6f2SLois Curfman McInnes     }
1034a9815358SBarry Smith   }
1035f3ef73ceSBarry Smith 
103629bbc08cSBarry Smith   if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
103729bbc08cSBarry Smith   if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first");
103829bbc08cSBarry Smith   if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option");
1039a8c6a408SBarry Smith   if (snes->vec_func == snes->vec_sol) {
104029bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector");
1041a8c6a408SBarry Smith   }
1042a703fe33SLois Curfman McInnes 
1043a703fe33SLois Curfman McInnes   /* Set the KSP stopping criterion to use the Eisenstat-Walker method */
10444b27c08aSLois Curfman McInnes   ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr);
10456831982aSBarry Smith   if (snes->ksp_ewconv && !iseqtr) {
104694b7f48cSBarry Smith     KSP ksp;
104794b7f48cSBarry Smith     ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
1048186905e3SBarry Smith     ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr);
1049a703fe33SLois Curfman McInnes   }
10504b27c08aSLois Curfman McInnes 
1051a703fe33SLois Curfman McInnes   if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);}
105282bf6240SBarry Smith   snes->setupcalled = 1;
10533a40ed3dSBarry Smith   PetscFunctionReturn(0);
10549b94acceSBarry Smith }
10559b94acceSBarry Smith 
10564a2ae208SSatish Balay #undef __FUNCT__
10574a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy"
10589b94acceSBarry Smith /*@C
10599b94acceSBarry Smith    SNESDestroy - Destroys the nonlinear solver context that was created
10609b94acceSBarry Smith    with SNESCreate().
10619b94acceSBarry Smith 
1062c7afd0dbSLois Curfman McInnes    Collective on SNES
1063c7afd0dbSLois Curfman McInnes 
10649b94acceSBarry Smith    Input Parameter:
10659b94acceSBarry Smith .  snes - the SNES context
10669b94acceSBarry Smith 
106736851e7fSLois Curfman McInnes    Level: beginner
106836851e7fSLois Curfman McInnes 
10699b94acceSBarry Smith .keywords: SNES, nonlinear, destroy
10709b94acceSBarry Smith 
107163a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve()
10729b94acceSBarry Smith @*/
1073dfbe8321SBarry Smith PetscErrorCode SNESDestroy(SNES snes)
10749b94acceSBarry Smith {
1075*6849ba73SBarry Smith   int i;
1076*6849ba73SBarry Smith   PetscErrorCode ierr;
10773a40ed3dSBarry Smith 
10783a40ed3dSBarry Smith   PetscFunctionBegin;
10794482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
10803a40ed3dSBarry Smith   if (--snes->refct > 0) PetscFunctionReturn(0);
1081d4bb536fSBarry Smith 
1082be0abb6dSBarry Smith   /* if memory was published with AMS then destroy it */
10830f5bd95cSBarry Smith   ierr = PetscObjectDepublish(snes);CHKERRQ(ierr);
1084be0abb6dSBarry Smith 
1085e1311b90SBarry Smith   if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);}
1086606d414cSSatish Balay   if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);}
10873a7fca6bSBarry Smith   if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);}
10883a7fca6bSBarry Smith   if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);}
10893ab0aad5SBarry Smith   if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);}
109094b7f48cSBarry Smith   ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);
1091522c5e43SBarry Smith   if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);}
1092b8a78c4aSBarry Smith   for (i=0; i<snes->numbermonitors; i++) {
1093b8a78c4aSBarry Smith     if (snes->monitordestroy[i]) {
1094b8a78c4aSBarry Smith       ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr);
1095b8a78c4aSBarry Smith     }
1096b8a78c4aSBarry Smith   }
1097b0a32e0cSBarry Smith   PetscLogObjectDestroy((PetscObject)snes);
10980452661fSBarry Smith   PetscHeaderDestroy((PetscObject)snes);
10993a40ed3dSBarry Smith   PetscFunctionReturn(0);
11009b94acceSBarry Smith }
11019b94acceSBarry Smith 
11029b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */
11039b94acceSBarry Smith 
11044a2ae208SSatish Balay #undef __FUNCT__
11054a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances"
11069b94acceSBarry Smith /*@
1107d7a720efSLois Curfman McInnes    SNESSetTolerances - Sets various parameters used in convergence tests.
11089b94acceSBarry Smith 
1109c7afd0dbSLois Curfman McInnes    Collective on SNES
1110c7afd0dbSLois Curfman McInnes 
11119b94acceSBarry Smith    Input Parameters:
1112c7afd0dbSLois Curfman McInnes +  snes - the SNES context
111333174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
111433174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
111533174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
111633174efeSLois Curfman McInnes            of the change in the solution between steps
111733174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1118c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1119fee21e36SBarry Smith 
112033174efeSLois Curfman McInnes    Options Database Keys:
1121c7afd0dbSLois Curfman McInnes +    -snes_atol <atol> - Sets atol
1122c7afd0dbSLois Curfman McInnes .    -snes_rtol <rtol> - Sets rtol
1123c7afd0dbSLois Curfman McInnes .    -snes_stol <stol> - Sets stol
1124c7afd0dbSLois Curfman McInnes .    -snes_max_it <maxit> - Sets maxit
1125c7afd0dbSLois Curfman McInnes -    -snes_max_funcs <maxf> - Sets maxf
11269b94acceSBarry Smith 
1127d7a720efSLois Curfman McInnes    Notes:
11289b94acceSBarry Smith    The default maximum number of iterations is 50.
11299b94acceSBarry Smith    The default maximum number of function evaluations is 1000.
11309b94acceSBarry Smith 
113136851e7fSLois Curfman McInnes    Level: intermediate
113236851e7fSLois Curfman McInnes 
113333174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances
11349b94acceSBarry Smith 
1135d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance()
11369b94acceSBarry Smith @*/
1137dfbe8321SBarry Smith PetscErrorCode SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf)
11389b94acceSBarry Smith {
11393a40ed3dSBarry Smith   PetscFunctionBegin;
11404482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1141d7a720efSLois Curfman McInnes   if (atol != PETSC_DEFAULT)  snes->atol      = atol;
1142d7a720efSLois Curfman McInnes   if (rtol != PETSC_DEFAULT)  snes->rtol      = rtol;
1143d7a720efSLois Curfman McInnes   if (stol != PETSC_DEFAULT)  snes->xtol      = stol;
1144d7a720efSLois Curfman McInnes   if (maxit != PETSC_DEFAULT) snes->max_its   = maxit;
1145d7a720efSLois Curfman McInnes   if (maxf != PETSC_DEFAULT)  snes->max_funcs = maxf;
11463a40ed3dSBarry Smith   PetscFunctionReturn(0);
11479b94acceSBarry Smith }
11489b94acceSBarry Smith 
11494a2ae208SSatish Balay #undef __FUNCT__
11504a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances"
11519b94acceSBarry Smith /*@
115233174efeSLois Curfman McInnes    SNESGetTolerances - Gets various parameters used in convergence tests.
115333174efeSLois Curfman McInnes 
1154c7afd0dbSLois Curfman McInnes    Not Collective
1155c7afd0dbSLois Curfman McInnes 
115633174efeSLois Curfman McInnes    Input Parameters:
1157c7afd0dbSLois Curfman McInnes +  snes - the SNES context
115833174efeSLois Curfman McInnes .  atol - absolute convergence tolerance
115933174efeSLois Curfman McInnes .  rtol - relative convergence tolerance
116033174efeSLois Curfman McInnes .  stol -  convergence tolerance in terms of the norm
116133174efeSLois Curfman McInnes            of the change in the solution between steps
116233174efeSLois Curfman McInnes .  maxit - maximum number of iterations
1163c7afd0dbSLois Curfman McInnes -  maxf - maximum number of function evaluations
1164fee21e36SBarry Smith 
116533174efeSLois Curfman McInnes    Notes:
116633174efeSLois Curfman McInnes    The user can specify PETSC_NULL for any parameter that is not needed.
116733174efeSLois Curfman McInnes 
116836851e7fSLois Curfman McInnes    Level: intermediate
116936851e7fSLois Curfman McInnes 
117033174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances
117133174efeSLois Curfman McInnes 
117233174efeSLois Curfman McInnes .seealso: SNESSetTolerances()
117333174efeSLois Curfman McInnes @*/
1174dfbe8321SBarry Smith PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf)
117533174efeSLois Curfman McInnes {
11763a40ed3dSBarry Smith   PetscFunctionBegin;
11774482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
117833174efeSLois Curfman McInnes   if (atol)  *atol  = snes->atol;
117933174efeSLois Curfman McInnes   if (rtol)  *rtol  = snes->rtol;
118033174efeSLois Curfman McInnes   if (stol)  *stol  = snes->xtol;
118133174efeSLois Curfman McInnes   if (maxit) *maxit = snes->max_its;
118233174efeSLois Curfman McInnes   if (maxf)  *maxf  = snes->max_funcs;
11833a40ed3dSBarry Smith   PetscFunctionReturn(0);
118433174efeSLois Curfman McInnes }
118533174efeSLois Curfman McInnes 
11864a2ae208SSatish Balay #undef __FUNCT__
11874a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance"
118833174efeSLois Curfman McInnes /*@
11899b94acceSBarry Smith    SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance.
11909b94acceSBarry Smith 
1191fee21e36SBarry Smith    Collective on SNES
1192fee21e36SBarry Smith 
1193c7afd0dbSLois Curfman McInnes    Input Parameters:
1194c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1195c7afd0dbSLois Curfman McInnes -  tol - tolerance
1196c7afd0dbSLois Curfman McInnes 
11979b94acceSBarry Smith    Options Database Key:
1198c7afd0dbSLois Curfman McInnes .  -snes_trtol <tol> - Sets tol
11999b94acceSBarry Smith 
120036851e7fSLois Curfman McInnes    Level: intermediate
120136851e7fSLois Curfman McInnes 
12029b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance
12039b94acceSBarry Smith 
1204d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance()
12059b94acceSBarry Smith @*/
1206dfbe8321SBarry Smith PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol)
12079b94acceSBarry Smith {
12083a40ed3dSBarry Smith   PetscFunctionBegin;
12094482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
12109b94acceSBarry Smith   snes->deltatol = tol;
12113a40ed3dSBarry Smith   PetscFunctionReturn(0);
12129b94acceSBarry Smith }
12139b94acceSBarry Smith 
1214df9fa365SBarry Smith /*
1215df9fa365SBarry Smith    Duplicate the lg monitors for SNES from KSP; for some reason with
1216df9fa365SBarry Smith    dynamic libraries things don't work under Sun4 if we just use
1217df9fa365SBarry Smith    macros instead of functions
1218df9fa365SBarry Smith */
12194a2ae208SSatish Balay #undef __FUNCT__
12204a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor"
1221dfbe8321SBarry Smith PetscErrorCode SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx)
1222ce1608b8SBarry Smith {
1223dfbe8321SBarry Smith   PetscErrorCode ierr;
1224ce1608b8SBarry Smith 
1225ce1608b8SBarry Smith   PetscFunctionBegin;
12264482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1227ce1608b8SBarry Smith   ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr);
1228ce1608b8SBarry Smith   PetscFunctionReturn(0);
1229ce1608b8SBarry Smith }
1230ce1608b8SBarry Smith 
12314a2ae208SSatish Balay #undef __FUNCT__
12324a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate"
1233dfbe8321SBarry Smith PetscErrorCode SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw)
1234df9fa365SBarry Smith {
1235dfbe8321SBarry Smith   PetscErrorCode ierr;
1236df9fa365SBarry Smith 
1237df9fa365SBarry Smith   PetscFunctionBegin;
1238df9fa365SBarry Smith   ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr);
1239df9fa365SBarry Smith   PetscFunctionReturn(0);
1240df9fa365SBarry Smith }
1241df9fa365SBarry Smith 
12424a2ae208SSatish Balay #undef __FUNCT__
12434a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy"
1244dfbe8321SBarry Smith PetscErrorCode SNESLGMonitorDestroy(PetscDrawLG draw)
1245df9fa365SBarry Smith {
1246dfbe8321SBarry Smith   PetscErrorCode ierr;
1247df9fa365SBarry Smith 
1248df9fa365SBarry Smith   PetscFunctionBegin;
1249df9fa365SBarry Smith   ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr);
1250df9fa365SBarry Smith   PetscFunctionReturn(0);
1251df9fa365SBarry Smith }
1252df9fa365SBarry Smith 
12539b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */
12549b94acceSBarry Smith 
12554a2ae208SSatish Balay #undef __FUNCT__
12564a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor"
12579b94acceSBarry Smith /*@C
12585cd90555SBarry Smith    SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every
12599b94acceSBarry Smith    iteration of the nonlinear solver to display the iteration's
12609b94acceSBarry Smith    progress.
12619b94acceSBarry Smith 
1262fee21e36SBarry Smith    Collective on SNES
1263fee21e36SBarry Smith 
1264c7afd0dbSLois Curfman McInnes    Input Parameters:
1265c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1266c7afd0dbSLois Curfman McInnes .  func - monitoring routine
1267b8a78c4aSBarry Smith .  mctx - [optional] user-defined context for private data for the
1268b3006f0bSLois Curfman McInnes           monitor routine (use PETSC_NULL if no context is desitre)
1269b3006f0bSLois Curfman McInnes -  monitordestroy - [optional] routine that frees monitor context
1270b3006f0bSLois Curfman McInnes           (may be PETSC_NULL)
12719b94acceSBarry Smith 
1272c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1273329f5518SBarry Smith $     int func(SNES snes,int its, PetscReal norm,void *mctx)
1274c7afd0dbSLois Curfman McInnes 
1275c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1276c7afd0dbSLois Curfman McInnes .    its - iteration number
1277c7afd0dbSLois Curfman McInnes .    norm - 2-norm function value (may be estimated)
127840a0c1c6SLois Curfman McInnes -    mctx - [optional] monitoring context
12799b94acceSBarry Smith 
12809665c990SLois Curfman McInnes    Options Database Keys:
1281c7afd0dbSLois Curfman McInnes +    -snes_monitor        - sets SNESDefaultMonitor()
1282c7afd0dbSLois Curfman McInnes .    -snes_xmonitor       - sets line graph monitor,
1283c7afd0dbSLois Curfman McInnes                             uses SNESLGMonitorCreate()
1284c7afd0dbSLois Curfman McInnes _    -snes_cancelmonitors - cancels all monitors that have
1285c7afd0dbSLois Curfman McInnes                             been hardwired into a code by
1286c7afd0dbSLois Curfman McInnes                             calls to SNESSetMonitor(), but
1287c7afd0dbSLois Curfman McInnes                             does not cancel those set via
1288c7afd0dbSLois Curfman McInnes                             the options database.
12899665c990SLois Curfman McInnes 
1290639f9d9dSBarry Smith    Notes:
12916bc08f3fSLois Curfman McInnes    Several different monitoring routines may be set by calling
12926bc08f3fSLois Curfman McInnes    SNESSetMonitor() multiple times; all will be called in the
12936bc08f3fSLois Curfman McInnes    order in which they were set.
1294639f9d9dSBarry Smith 
129536851e7fSLois Curfman McInnes    Level: intermediate
129636851e7fSLois Curfman McInnes 
12979b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor
12989b94acceSBarry Smith 
12995cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor()
13009b94acceSBarry Smith @*/
1301*6849ba73SBarry Smith PetscErrorCode SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,int,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*))
13029b94acceSBarry Smith {
13033a40ed3dSBarry Smith   PetscFunctionBegin;
13044482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1305639f9d9dSBarry Smith   if (snes->numbermonitors >= MAXSNESMONITORS) {
130629bbc08cSBarry Smith     SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set");
1307639f9d9dSBarry Smith   }
1308639f9d9dSBarry Smith 
1309639f9d9dSBarry Smith   snes->monitor[snes->numbermonitors]           = func;
1310b8a78c4aSBarry Smith   snes->monitordestroy[snes->numbermonitors]    = monitordestroy;
1311639f9d9dSBarry Smith   snes->monitorcontext[snes->numbermonitors++]  = (void*)mctx;
13123a40ed3dSBarry Smith   PetscFunctionReturn(0);
13139b94acceSBarry Smith }
13149b94acceSBarry Smith 
13154a2ae208SSatish Balay #undef __FUNCT__
13164a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor"
13175cd90555SBarry Smith /*@C
13185cd90555SBarry Smith    SNESClearMonitor - Clears all the monitor functions for a SNES object.
13195cd90555SBarry Smith 
1320c7afd0dbSLois Curfman McInnes    Collective on SNES
1321c7afd0dbSLois Curfman McInnes 
13225cd90555SBarry Smith    Input Parameters:
13235cd90555SBarry Smith .  snes - the SNES context
13245cd90555SBarry Smith 
13251a480d89SAdministrator    Options Database Key:
1326c7afd0dbSLois Curfman McInnes .  -snes_cancelmonitors - cancels all monitors that have been hardwired
1327c7afd0dbSLois Curfman McInnes     into a code by calls to SNESSetMonitor(), but does not cancel those
1328c7afd0dbSLois Curfman McInnes     set via the options database
13295cd90555SBarry Smith 
13305cd90555SBarry Smith    Notes:
13315cd90555SBarry Smith    There is no way to clear one specific monitor from a SNES object.
13325cd90555SBarry Smith 
133336851e7fSLois Curfman McInnes    Level: intermediate
133436851e7fSLois Curfman McInnes 
13355cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor
13365cd90555SBarry Smith 
13375cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor()
13385cd90555SBarry Smith @*/
1339dfbe8321SBarry Smith PetscErrorCode SNESClearMonitor(SNES snes)
13405cd90555SBarry Smith {
13415cd90555SBarry Smith   PetscFunctionBegin;
13424482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
13435cd90555SBarry Smith   snes->numbermonitors = 0;
13445cd90555SBarry Smith   PetscFunctionReturn(0);
13455cd90555SBarry Smith }
13465cd90555SBarry Smith 
13474a2ae208SSatish Balay #undef __FUNCT__
13484a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest"
13499b94acceSBarry Smith /*@C
13509b94acceSBarry Smith    SNESSetConvergenceTest - Sets the function that is to be used
13519b94acceSBarry Smith    to test for convergence of the nonlinear iterative solution.
13529b94acceSBarry Smith 
1353fee21e36SBarry Smith    Collective on SNES
1354fee21e36SBarry Smith 
1355c7afd0dbSLois Curfman McInnes    Input Parameters:
1356c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1357c7afd0dbSLois Curfman McInnes .  func - routine to test for convergence
1358c7afd0dbSLois Curfman McInnes -  cctx - [optional] context for private data for the convergence routine
1359c7afd0dbSLois Curfman McInnes           (may be PETSC_NULL)
13609b94acceSBarry Smith 
1361c7afd0dbSLois Curfman McInnes    Calling sequence of func:
1362*6849ba73SBarry Smith $     PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx)
1363c7afd0dbSLois Curfman McInnes 
1364c7afd0dbSLois Curfman McInnes +    snes - the SNES context
1365c7afd0dbSLois Curfman McInnes .    cctx - [optional] convergence context
1366184914b5SBarry Smith .    reason - reason for convergence/divergence
1367c7afd0dbSLois Curfman McInnes .    xnorm - 2-norm of current iterate
13684b27c08aSLois Curfman McInnes .    gnorm - 2-norm of current step
13694b27c08aSLois Curfman McInnes -    f - 2-norm of function
13709b94acceSBarry Smith 
137136851e7fSLois Curfman McInnes    Level: advanced
137236851e7fSLois Curfman McInnes 
13739b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test
13749b94acceSBarry Smith 
13754b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR()
13769b94acceSBarry Smith @*/
1377*6849ba73SBarry Smith PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx)
13789b94acceSBarry Smith {
13793a40ed3dSBarry Smith   PetscFunctionBegin;
13804482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
13819b94acceSBarry Smith   (snes)->converged = func;
13829b94acceSBarry Smith   (snes)->cnvP      = cctx;
13833a40ed3dSBarry Smith   PetscFunctionReturn(0);
13849b94acceSBarry Smith }
13859b94acceSBarry Smith 
13864a2ae208SSatish Balay #undef __FUNCT__
13874a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason"
1388184914b5SBarry Smith /*@C
1389184914b5SBarry Smith    SNESGetConvergedReason - Gets the reason the SNES iteration was stopped.
1390184914b5SBarry Smith 
1391184914b5SBarry Smith    Not Collective
1392184914b5SBarry Smith 
1393184914b5SBarry Smith    Input Parameter:
1394184914b5SBarry Smith .  snes - the SNES context
1395184914b5SBarry Smith 
1396184914b5SBarry Smith    Output Parameter:
1397e090d566SSatish Balay .  reason - negative value indicates diverged, positive value converged, see petscsnes.h or the
1398184914b5SBarry Smith             manual pages for the individual convergence tests for complete lists
1399184914b5SBarry Smith 
1400184914b5SBarry Smith    Level: intermediate
1401184914b5SBarry Smith 
1402184914b5SBarry Smith    Notes: Can only be called after the call the SNESSolve() is complete.
1403184914b5SBarry Smith 
1404184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test
1405184914b5SBarry Smith 
14064b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason
1407184914b5SBarry Smith @*/
1408dfbe8321SBarry Smith PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason)
1409184914b5SBarry Smith {
1410184914b5SBarry Smith   PetscFunctionBegin;
14114482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14124482741eSBarry Smith   PetscValidPointer(reason,2);
1413184914b5SBarry Smith   *reason = snes->reason;
1414184914b5SBarry Smith   PetscFunctionReturn(0);
1415184914b5SBarry Smith }
1416184914b5SBarry Smith 
14174a2ae208SSatish Balay #undef __FUNCT__
14184a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory"
1419c9005455SLois Curfman McInnes /*@
1420c9005455SLois Curfman McInnes    SNESSetConvergenceHistory - Sets the array used to hold the convergence history.
1421c9005455SLois Curfman McInnes 
1422fee21e36SBarry Smith    Collective on SNES
1423fee21e36SBarry Smith 
1424c7afd0dbSLois Curfman McInnes    Input Parameters:
1425c7afd0dbSLois Curfman McInnes +  snes - iterative context obtained from SNESCreate()
1426c7afd0dbSLois Curfman McInnes .  a   - array to hold history
1427cd5578b5SBarry Smith .  its - integer array holds the number of linear iterations for each solve.
1428758f92a0SBarry Smith .  na  - size of a and its
142964731454SLois Curfman McInnes -  reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero,
1430758f92a0SBarry Smith            else it continues storing new values for new nonlinear solves after the old ones
1431c7afd0dbSLois Curfman McInnes 
1432c9005455SLois Curfman McInnes    Notes:
14334b27c08aSLois Curfman McInnes    If set, this array will contain the function norms computed
1434c9005455SLois Curfman McInnes    at each step.
1435c9005455SLois Curfman McInnes 
1436c9005455SLois Curfman McInnes    This routine is useful, e.g., when running a code for purposes
1437c9005455SLois Curfman McInnes    of accurate performance monitoring, when no I/O should be done
1438c9005455SLois Curfman McInnes    during the section of code that is being timed.
1439c9005455SLois Curfman McInnes 
144036851e7fSLois Curfman McInnes    Level: intermediate
144136851e7fSLois Curfman McInnes 
1442c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history
1443758f92a0SBarry Smith 
144408405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory()
1445758f92a0SBarry Smith 
1446c9005455SLois Curfman McInnes @*/
1447dfbe8321SBarry Smith PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],int *its,int na,PetscTruth reset)
1448c9005455SLois Curfman McInnes {
14493a40ed3dSBarry Smith   PetscFunctionBegin;
14504482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
14514482741eSBarry Smith   if (na) PetscValidScalarPointer(a,2);
1452c9005455SLois Curfman McInnes   snes->conv_hist       = a;
1453758f92a0SBarry Smith   snes->conv_hist_its   = its;
1454758f92a0SBarry Smith   snes->conv_hist_max   = na;
1455758f92a0SBarry Smith   snes->conv_hist_reset = reset;
1456758f92a0SBarry Smith   PetscFunctionReturn(0);
1457758f92a0SBarry Smith }
1458758f92a0SBarry Smith 
14594a2ae208SSatish Balay #undef __FUNCT__
14604a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory"
14610c4c9dddSBarry Smith /*@C
1462758f92a0SBarry Smith    SNESGetConvergenceHistory - Gets the array used to hold the convergence history.
1463758f92a0SBarry Smith 
1464758f92a0SBarry Smith    Collective on SNES
1465758f92a0SBarry Smith 
1466758f92a0SBarry Smith    Input Parameter:
1467758f92a0SBarry Smith .  snes - iterative context obtained from SNESCreate()
1468758f92a0SBarry Smith 
1469758f92a0SBarry Smith    Output Parameters:
1470758f92a0SBarry Smith .  a   - array to hold history
1471758f92a0SBarry Smith .  its - integer array holds the number of linear iterations (or
1472758f92a0SBarry Smith          negative if not converged) for each solve.
1473758f92a0SBarry Smith -  na  - size of a and its
1474758f92a0SBarry Smith 
1475758f92a0SBarry Smith    Notes:
1476758f92a0SBarry Smith     The calling sequence for this routine in Fortran is
1477758f92a0SBarry Smith $   call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr)
1478758f92a0SBarry Smith 
1479758f92a0SBarry Smith    This routine is useful, e.g., when running a code for purposes
1480758f92a0SBarry Smith    of accurate performance monitoring, when no I/O should be done
1481758f92a0SBarry Smith    during the section of code that is being timed.
1482758f92a0SBarry Smith 
1483758f92a0SBarry Smith    Level: intermediate
1484758f92a0SBarry Smith 
1485758f92a0SBarry Smith .keywords: SNES, get, convergence, history
1486758f92a0SBarry Smith 
1487758f92a0SBarry Smith .seealso: SNESSetConvergencHistory()
1488758f92a0SBarry Smith 
1489758f92a0SBarry Smith @*/
1490dfbe8321SBarry Smith PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],int *its[],int *na)
1491758f92a0SBarry Smith {
1492758f92a0SBarry Smith   PetscFunctionBegin;
14934482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
1494758f92a0SBarry Smith   if (a)   *a   = snes->conv_hist;
1495758f92a0SBarry Smith   if (its) *its = snes->conv_hist_its;
1496758f92a0SBarry Smith   if (na) *na   = snes->conv_hist_len;
14973a40ed3dSBarry Smith   PetscFunctionReturn(0);
1498c9005455SLois Curfman McInnes }
1499c9005455SLois Curfman McInnes 
1500e74ef692SMatthew Knepley #undef __FUNCT__
1501e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC"
1502ac226902SBarry Smith /*@C
150376b2cf59SMatthew Knepley   SNESSetRhsBC - Sets the function which applies boundary conditions
150476b2cf59SMatthew Knepley   to the Rhs of each system.
150576b2cf59SMatthew Knepley 
150676b2cf59SMatthew Knepley   Collective on SNES
150776b2cf59SMatthew Knepley 
150876b2cf59SMatthew Knepley   Input Parameters:
150976b2cf59SMatthew Knepley . snes - The nonlinear solver context
151076b2cf59SMatthew Knepley . func - The function
151176b2cf59SMatthew Knepley 
151276b2cf59SMatthew Knepley   Calling sequence of func:
151376b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx);
151476b2cf59SMatthew Knepley 
151576b2cf59SMatthew Knepley . rhs - The current rhs vector
151676b2cf59SMatthew Knepley . ctx - The user-context
151776b2cf59SMatthew Knepley 
151876b2cf59SMatthew Knepley   Level: intermediate
151976b2cf59SMatthew Knepley 
152076b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions
152176b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate()
152276b2cf59SMatthew Knepley @*/
1523*6849ba73SBarry Smith PetscErrorCode SNESSetRhsBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *))
152476b2cf59SMatthew Knepley {
152576b2cf59SMatthew Knepley   PetscFunctionBegin;
15264482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
152776b2cf59SMatthew Knepley   snes->applyrhsbc = func;
152876b2cf59SMatthew Knepley   PetscFunctionReturn(0);
152976b2cf59SMatthew Knepley }
153076b2cf59SMatthew Knepley 
1531e74ef692SMatthew Knepley #undef __FUNCT__
1532e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC"
153376b2cf59SMatthew Knepley /*@
153476b2cf59SMatthew Knepley   SNESDefaultRhsBC - The default boundary condition function which does nothing.
153576b2cf59SMatthew Knepley 
153676b2cf59SMatthew Knepley   Not collective
153776b2cf59SMatthew Knepley 
153876b2cf59SMatthew Knepley   Input Parameters:
153976b2cf59SMatthew Knepley . snes - The nonlinear solver context
154076b2cf59SMatthew Knepley . rhs  - The Rhs
154176b2cf59SMatthew Knepley . ctx  - The user-context
154276b2cf59SMatthew Knepley 
154376b2cf59SMatthew Knepley   Level: beginner
154476b2cf59SMatthew Knepley 
154576b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions
154676b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate()
154776b2cf59SMatthew Knepley @*/
1548dfbe8321SBarry Smith PetscErrorCode SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx)
154976b2cf59SMatthew Knepley {
155076b2cf59SMatthew Knepley   PetscFunctionBegin;
155176b2cf59SMatthew Knepley   PetscFunctionReturn(0);
155276b2cf59SMatthew Knepley }
155376b2cf59SMatthew Knepley 
1554e74ef692SMatthew Knepley #undef __FUNCT__
1555e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC"
1556ac226902SBarry Smith /*@C
155776b2cf59SMatthew Knepley   SNESSetSolutionBC - Sets the function which applies boundary conditions
155876b2cf59SMatthew Knepley   to the solution of each system.
155976b2cf59SMatthew Knepley 
156076b2cf59SMatthew Knepley   Collective on SNES
156176b2cf59SMatthew Knepley 
156276b2cf59SMatthew Knepley   Input Parameters:
156376b2cf59SMatthew Knepley . snes - The nonlinear solver context
156476b2cf59SMatthew Knepley . func - The function
156576b2cf59SMatthew Knepley 
156676b2cf59SMatthew Knepley   Calling sequence of func:
15679d06fb6bSMatthew Knepley . func (SNES snes, Vec rsol, void *ctx);
156876b2cf59SMatthew Knepley 
156976b2cf59SMatthew Knepley . sol - The current solution vector
157076b2cf59SMatthew Knepley . ctx - The user-context
157176b2cf59SMatthew Knepley 
157276b2cf59SMatthew Knepley   Level: intermediate
157376b2cf59SMatthew Knepley 
157476b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions
157576b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate()
157676b2cf59SMatthew Knepley @*/
1577*6849ba73SBarry Smith PetscErrorCode SNESSetSolutionBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *))
157876b2cf59SMatthew Knepley {
157976b2cf59SMatthew Knepley   PetscFunctionBegin;
15804482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
158176b2cf59SMatthew Knepley   snes->applysolbc = func;
158276b2cf59SMatthew Knepley   PetscFunctionReturn(0);
158376b2cf59SMatthew Knepley }
158476b2cf59SMatthew Knepley 
1585e74ef692SMatthew Knepley #undef __FUNCT__
1586e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC"
158776b2cf59SMatthew Knepley /*@
158876b2cf59SMatthew Knepley   SNESDefaultSolutionBC - The default boundary condition function which does nothing.
158976b2cf59SMatthew Knepley 
159076b2cf59SMatthew Knepley   Not collective
159176b2cf59SMatthew Knepley 
159276b2cf59SMatthew Knepley   Input Parameters:
159376b2cf59SMatthew Knepley . snes - The nonlinear solver context
159476b2cf59SMatthew Knepley . sol  - The solution
159576b2cf59SMatthew Knepley . ctx  - The user-context
159676b2cf59SMatthew Knepley 
159776b2cf59SMatthew Knepley   Level: beginner
159876b2cf59SMatthew Knepley 
159976b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions
160076b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate()
160176b2cf59SMatthew Knepley @*/
1602dfbe8321SBarry Smith PetscErrorCode SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx)
160376b2cf59SMatthew Knepley {
160476b2cf59SMatthew Knepley   PetscFunctionBegin;
160576b2cf59SMatthew Knepley   PetscFunctionReturn(0);
160676b2cf59SMatthew Knepley }
160776b2cf59SMatthew Knepley 
1608e74ef692SMatthew Knepley #undef __FUNCT__
1609e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate"
1610ac226902SBarry Smith /*@C
161176b2cf59SMatthew Knepley   SNESSetUpdate - Sets the general-purpose update function called
161276b2cf59SMatthew Knepley   at the beginning of every step of the iteration.
161376b2cf59SMatthew Knepley 
161476b2cf59SMatthew Knepley   Collective on SNES
161576b2cf59SMatthew Knepley 
161676b2cf59SMatthew Knepley   Input Parameters:
161776b2cf59SMatthew Knepley . snes - The nonlinear solver context
161876b2cf59SMatthew Knepley . func - The function
161976b2cf59SMatthew Knepley 
162076b2cf59SMatthew Knepley   Calling sequence of func:
162176b2cf59SMatthew Knepley . func (TS ts, int step);
162276b2cf59SMatthew Knepley 
162376b2cf59SMatthew Knepley . step - The current step of the iteration
162476b2cf59SMatthew Knepley 
162576b2cf59SMatthew Knepley   Level: intermediate
162676b2cf59SMatthew Knepley 
162776b2cf59SMatthew Knepley .keywords: SNES, update
162876b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC()
162976b2cf59SMatthew Knepley @*/
1630*6849ba73SBarry Smith PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, int))
163176b2cf59SMatthew Knepley {
163276b2cf59SMatthew Knepley   PetscFunctionBegin;
16334482741eSBarry Smith   PetscValidHeaderSpecific(snes, SNES_COOKIE,1);
163476b2cf59SMatthew Knepley   snes->update = func;
163576b2cf59SMatthew Knepley   PetscFunctionReturn(0);
163676b2cf59SMatthew Knepley }
163776b2cf59SMatthew Knepley 
1638e74ef692SMatthew Knepley #undef __FUNCT__
1639e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate"
164076b2cf59SMatthew Knepley /*@
164176b2cf59SMatthew Knepley   SNESDefaultUpdate - The default update function which does nothing.
164276b2cf59SMatthew Knepley 
164376b2cf59SMatthew Knepley   Not collective
164476b2cf59SMatthew Knepley 
164576b2cf59SMatthew Knepley   Input Parameters:
164676b2cf59SMatthew Knepley . snes - The nonlinear solver context
164776b2cf59SMatthew Knepley . step - The current step of the iteration
164876b2cf59SMatthew Knepley 
1649205452f4SMatthew Knepley   Level: intermediate
1650205452f4SMatthew Knepley 
165176b2cf59SMatthew Knepley .keywords: SNES, update
165276b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC()
165376b2cf59SMatthew Knepley @*/
1654dfbe8321SBarry Smith PetscErrorCode SNESDefaultUpdate(SNES snes, int step)
165576b2cf59SMatthew Knepley {
165676b2cf59SMatthew Knepley   PetscFunctionBegin;
165776b2cf59SMatthew Knepley   PetscFunctionReturn(0);
165876b2cf59SMatthew Knepley }
165976b2cf59SMatthew Knepley 
16604a2ae208SSatish Balay #undef __FUNCT__
16614a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private"
16629b94acceSBarry Smith /*
16639b94acceSBarry Smith    SNESScaleStep_Private - Scales a step so that its length is less than the
16649b94acceSBarry Smith    positive parameter delta.
16659b94acceSBarry Smith 
16669b94acceSBarry Smith     Input Parameters:
1667c7afd0dbSLois Curfman McInnes +   snes - the SNES context
16689b94acceSBarry Smith .   y - approximate solution of linear system
16699b94acceSBarry Smith .   fnorm - 2-norm of current function
1670c7afd0dbSLois Curfman McInnes -   delta - trust region size
16719b94acceSBarry Smith 
16729b94acceSBarry Smith     Output Parameters:
1673c7afd0dbSLois Curfman McInnes +   gpnorm - predicted function norm at the new point, assuming local
16749b94acceSBarry Smith     linearization.  The value is zero if the step lies within the trust
16759b94acceSBarry Smith     region, and exceeds zero otherwise.
1676c7afd0dbSLois Curfman McInnes -   ynorm - 2-norm of the step
16779b94acceSBarry Smith 
16789b94acceSBarry Smith     Note:
16794b27c08aSLois Curfman McInnes     For non-trust region methods such as SNESLS, the parameter delta
16809b94acceSBarry Smith     is set to be the maximum allowable step size.
16819b94acceSBarry Smith 
16829b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step
16839b94acceSBarry Smith */
1684dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm)
16859b94acceSBarry Smith {
1686064f8208SBarry Smith   PetscReal   nrm;
1687ea709b57SSatish Balay   PetscScalar cnorm;
1688dfbe8321SBarry Smith   PetscErrorCode ierr;
16893a40ed3dSBarry Smith 
16903a40ed3dSBarry Smith   PetscFunctionBegin;
16914482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
16924482741eSBarry Smith   PetscValidHeaderSpecific(y,VEC_COOKIE,2);
1693c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,y,2);
1694184914b5SBarry Smith 
1695064f8208SBarry Smith   ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr);
1696064f8208SBarry Smith   if (nrm > *delta) {
1697064f8208SBarry Smith      nrm = *delta/nrm;
1698064f8208SBarry Smith      *gpnorm = (1.0 - nrm)*(*fnorm);
1699064f8208SBarry Smith      cnorm = nrm;
1700329f5518SBarry Smith      ierr = VecScale(&cnorm,y);CHKERRQ(ierr);
17019b94acceSBarry Smith      *ynorm = *delta;
17029b94acceSBarry Smith   } else {
17039b94acceSBarry Smith      *gpnorm = 0.0;
1704064f8208SBarry Smith      *ynorm = nrm;
17059b94acceSBarry Smith   }
17063a40ed3dSBarry Smith   PetscFunctionReturn(0);
17079b94acceSBarry Smith }
17089b94acceSBarry Smith 
17095968eb51SBarry Smith static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero",
17105968eb51SBarry Smith                                          "not currently used",
17115968eb51SBarry Smith                                          "line search failed",
17125968eb51SBarry Smith                                          "reach maximum number of iterations",
17135968eb51SBarry Smith                                          "function norm became NaN (not a number)",
17145968eb51SBarry Smith                                          "not currently used",
17155968eb51SBarry Smith                                          "number of function computations exceeded",
17165968eb51SBarry Smith                                          "not currently used",
17175968eb51SBarry Smith                                          "still iterating",
17185968eb51SBarry Smith                                          "not currently used",
17195968eb51SBarry Smith                                          "absolute size of function norm",
17205968eb51SBarry Smith                                          "relative decrease in function norm",
17215968eb51SBarry Smith                                          "step size is small",
17225968eb51SBarry Smith                                          "not currently used",
17235968eb51SBarry Smith                                          "not currently used",
17245968eb51SBarry Smith                                          "small trust region"};
17255968eb51SBarry Smith 
17264a2ae208SSatish Balay #undef __FUNCT__
17274a2ae208SSatish Balay #define __FUNCT__ "SNESSolve"
17289b94acceSBarry Smith /*@
17299b94acceSBarry Smith    SNESSolve - Solves a nonlinear system.  Call SNESSolve after calling
173063a78c88SLois Curfman McInnes    SNESCreate() and optional routines of the form SNESSetXXX().
17319b94acceSBarry Smith 
1732c7afd0dbSLois Curfman McInnes    Collective on SNES
1733c7afd0dbSLois Curfman McInnes 
1734b2002411SLois Curfman McInnes    Input Parameters:
1735c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1736c7afd0dbSLois Curfman McInnes -  x - the solution vector
17379b94acceSBarry Smith 
1738b2002411SLois Curfman McInnes    Notes:
17398ddd3da0SLois Curfman McInnes    The user should initialize the vector,x, with the initial guess
17408ddd3da0SLois Curfman McInnes    for the nonlinear solve prior to calling SNESSolve.  In particular,
17418ddd3da0SLois Curfman McInnes    to employ an initial guess of zero, the user should explicitly set
17428ddd3da0SLois Curfman McInnes    this vector to zero by calling VecSet().
17438ddd3da0SLois Curfman McInnes 
174436851e7fSLois Curfman McInnes    Level: beginner
174536851e7fSLois Curfman McInnes 
17469b94acceSBarry Smith .keywords: SNES, nonlinear, solve
17479b94acceSBarry Smith 
17483ab0aad5SBarry Smith .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs()
17499b94acceSBarry Smith @*/
1750dfbe8321SBarry Smith PetscErrorCode SNESSolve(SNES snes,Vec x)
17519b94acceSBarry Smith {
1752dfbe8321SBarry Smith   PetscErrorCode ierr;
1753f1af5d2fSBarry Smith   PetscTruth flg;
1754052efed2SBarry Smith 
17553a40ed3dSBarry Smith   PetscFunctionBegin;
17564482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
17574482741eSBarry Smith   PetscValidHeaderSpecific(x,VEC_COOKIE,2);
1758c9780b6fSBarry Smith   PetscCheckSameComm(snes,1,x,2);
175929bbc08cSBarry Smith   if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()");
176074637425SBarry Smith 
176182bf6240SBarry Smith   if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);}
1762c4fc05e7SBarry Smith   else {snes->vec_sol = snes->vec_sol_always = x;}
1763758f92a0SBarry Smith   if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0;
1764d5ba7fb7SMatthew Knepley   ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
176550ffb88aSMatthew Knepley   snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0;
1766c9780b6fSBarry Smith   ierr = (*(snes)->solve)(snes);CHKERRQ(ierr);
1767d5ba7fb7SMatthew Knepley   ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr);
1768b0a32e0cSBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr);
17698b179ff8SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); }
1770da9b6338SBarry Smith   ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr);
1771da9b6338SBarry Smith   if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); }
17725968eb51SBarry Smith   if (snes->printreason) {
17735968eb51SBarry Smith     if (snes->reason > 0) {
17745968eb51SBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr);
17755968eb51SBarry Smith     } else {
17765968eb51SBarry Smith       ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr);
17775968eb51SBarry Smith     }
17785968eb51SBarry Smith   }
17795968eb51SBarry Smith 
17803a40ed3dSBarry Smith   PetscFunctionReturn(0);
17819b94acceSBarry Smith }
17829b94acceSBarry Smith 
17839b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */
17849b94acceSBarry Smith 
17854a2ae208SSatish Balay #undef __FUNCT__
17864a2ae208SSatish Balay #define __FUNCT__ "SNESSetType"
178782bf6240SBarry Smith /*@C
17884b0e389bSBarry Smith    SNESSetType - Sets the method for the nonlinear solver.
17899b94acceSBarry Smith 
1790fee21e36SBarry Smith    Collective on SNES
1791fee21e36SBarry Smith 
1792c7afd0dbSLois Curfman McInnes    Input Parameters:
1793c7afd0dbSLois Curfman McInnes +  snes - the SNES context
1794454a90a3SBarry Smith -  type - a known method
1795c7afd0dbSLois Curfman McInnes 
1796c7afd0dbSLois Curfman McInnes    Options Database Key:
1797454a90a3SBarry Smith .  -snes_type <type> - Sets the method; use -help for a list
1798c7afd0dbSLois Curfman McInnes    of available methods (for instance, ls or tr)
1799ae12b187SLois Curfman McInnes 
18009b94acceSBarry Smith    Notes:
1801e090d566SSatish Balay    See "petsc/include/petscsnes.h" for available methods (for instance)
18024b27c08aSLois Curfman McInnes +    SNESLS - Newton's method with line search
1803c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
18044b27c08aSLois Curfman McInnes .    SNESTR - Newton's method with trust region
1805c7afd0dbSLois Curfman McInnes      (systems of nonlinear equations)
18069b94acceSBarry Smith 
1807ae12b187SLois Curfman McInnes   Normally, it is best to use the SNESSetFromOptions() command and then
1808ae12b187SLois Curfman McInnes   set the SNES solver type from the options database rather than by using
1809ae12b187SLois Curfman McInnes   this routine.  Using the options database provides the user with
1810ae12b187SLois Curfman McInnes   maximum flexibility in evaluating the many nonlinear solvers.
1811ae12b187SLois Curfman McInnes   The SNESSetType() routine is provided for those situations where it
1812ae12b187SLois Curfman McInnes   is necessary to set the nonlinear solver independently of the command
1813ae12b187SLois Curfman McInnes   line or options database.  This might be the case, for example, when
1814ae12b187SLois Curfman McInnes   the choice of solver changes during the execution of the program,
1815ae12b187SLois Curfman McInnes   and the user's application is taking responsibility for choosing the
1816b0a32e0cSBarry Smith   appropriate method.
181736851e7fSLois Curfman McInnes 
181836851e7fSLois Curfman McInnes   Level: intermediate
1819a703fe33SLois Curfman McInnes 
1820454a90a3SBarry Smith .keywords: SNES, set, type
1821435da068SBarry Smith 
1822435da068SBarry Smith .seealso: SNESType, SNESCreate()
1823435da068SBarry Smith 
18249b94acceSBarry Smith @*/
1825dfbe8321SBarry Smith PetscErrorCode SNESSetType(SNES snes,const SNESType type)
18269b94acceSBarry Smith {
1827dfbe8321SBarry Smith   PetscErrorCode ierr,(*r)(SNES);
18286831982aSBarry Smith   PetscTruth match;
18293a40ed3dSBarry Smith 
18303a40ed3dSBarry Smith   PetscFunctionBegin;
18314482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
18324482741eSBarry Smith   PetscValidCharPointer(type,2);
183382bf6240SBarry Smith 
18346831982aSBarry Smith   ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr);
18350f5bd95cSBarry Smith   if (match) PetscFunctionReturn(0);
183682bf6240SBarry Smith 
183782bf6240SBarry Smith   if (snes->setupcalled) {
1838e1311b90SBarry Smith     ierr       = (*(snes)->destroy)(snes);CHKERRQ(ierr);
183982bf6240SBarry Smith     snes->data = 0;
184082bf6240SBarry Smith   }
18417d1a2b2bSBarry Smith 
18429b94acceSBarry Smith   /* Get the function pointers for the iterative method requested */
184382bf6240SBarry Smith   if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);}
1844b9617806SBarry Smith   ierr =  PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr);
1845958c9bccSBarry Smith   if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type);
1846606d414cSSatish Balay   if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);}
184782bf6240SBarry Smith   snes->data = 0;
18483a40ed3dSBarry Smith   ierr = (*r)(snes);CHKERRQ(ierr);
1849454a90a3SBarry Smith   ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr);
18503a40ed3dSBarry Smith   PetscFunctionReturn(0);
18519b94acceSBarry Smith }
18529b94acceSBarry Smith 
1853a847f771SSatish Balay 
18549b94acceSBarry Smith /* --------------------------------------------------------------------- */
18554a2ae208SSatish Balay #undef __FUNCT__
18564a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy"
18579b94acceSBarry Smith /*@C
18589b94acceSBarry Smith    SNESRegisterDestroy - Frees the list of nonlinear solvers that were
1859f1af5d2fSBarry Smith    registered by SNESRegisterDynamic().
18609b94acceSBarry Smith 
1861fee21e36SBarry Smith    Not Collective
1862fee21e36SBarry Smith 
186336851e7fSLois Curfman McInnes    Level: advanced
186436851e7fSLois Curfman McInnes 
18659b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy
18669b94acceSBarry Smith 
18679b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll()
18689b94acceSBarry Smith @*/
1869dfbe8321SBarry Smith PetscErrorCode SNESRegisterDestroy(void)
18709b94acceSBarry Smith {
1871dfbe8321SBarry Smith   PetscErrorCode ierr;
187282bf6240SBarry Smith 
18733a40ed3dSBarry Smith   PetscFunctionBegin;
187482bf6240SBarry Smith   if (SNESList) {
1875b0a32e0cSBarry Smith     ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr);
187682bf6240SBarry Smith     SNESList = 0;
18779b94acceSBarry Smith   }
18784c49b128SBarry Smith   SNESRegisterAllCalled = PETSC_FALSE;
18793a40ed3dSBarry Smith   PetscFunctionReturn(0);
18809b94acceSBarry Smith }
18819b94acceSBarry Smith 
18824a2ae208SSatish Balay #undef __FUNCT__
18834a2ae208SSatish Balay #define __FUNCT__ "SNESGetType"
18849b94acceSBarry Smith /*@C
18859a28b0a6SLois Curfman McInnes    SNESGetType - Gets the SNES method type and name (as a string).
18869b94acceSBarry Smith 
1887c7afd0dbSLois Curfman McInnes    Not Collective
1888c7afd0dbSLois Curfman McInnes 
18899b94acceSBarry Smith    Input Parameter:
18904b0e389bSBarry Smith .  snes - nonlinear solver context
18919b94acceSBarry Smith 
18929b94acceSBarry Smith    Output Parameter:
18933a7fca6bSBarry Smith .  type - SNES method (a character string)
18949b94acceSBarry Smith 
189536851e7fSLois Curfman McInnes    Level: intermediate
189636851e7fSLois Curfman McInnes 
1897454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name
18989b94acceSBarry Smith @*/
1899dfbe8321SBarry Smith PetscErrorCode SNESGetType(SNES snes,SNESType *type)
19009b94acceSBarry Smith {
19013a40ed3dSBarry Smith   PetscFunctionBegin;
19024482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19034482741eSBarry Smith   PetscValidPointer(type,2);
1904454a90a3SBarry Smith   *type = snes->type_name;
19053a40ed3dSBarry Smith   PetscFunctionReturn(0);
19069b94acceSBarry Smith }
19079b94acceSBarry Smith 
19084a2ae208SSatish Balay #undef __FUNCT__
19094a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution"
19109b94acceSBarry Smith /*@C
19119b94acceSBarry Smith    SNESGetSolution - Returns the vector where the approximate solution is
19129b94acceSBarry Smith    stored.
19139b94acceSBarry Smith 
1914c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1915c7afd0dbSLois Curfman McInnes 
19169b94acceSBarry Smith    Input Parameter:
19179b94acceSBarry Smith .  snes - the SNES context
19189b94acceSBarry Smith 
19199b94acceSBarry Smith    Output Parameter:
19209b94acceSBarry Smith .  x - the solution
19219b94acceSBarry Smith 
192236851e7fSLois Curfman McInnes    Level: advanced
192336851e7fSLois Curfman McInnes 
19249b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution
19259b94acceSBarry Smith 
19264b27c08aSLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetSolutionUpdate()
19279b94acceSBarry Smith @*/
1928dfbe8321SBarry Smith PetscErrorCode SNESGetSolution(SNES snes,Vec *x)
19299b94acceSBarry Smith {
19303a40ed3dSBarry Smith   PetscFunctionBegin;
19314482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19324482741eSBarry Smith   PetscValidPointer(x,2);
19339b94acceSBarry Smith   *x = snes->vec_sol_always;
19343a40ed3dSBarry Smith   PetscFunctionReturn(0);
19359b94acceSBarry Smith }
19369b94acceSBarry Smith 
19374a2ae208SSatish Balay #undef __FUNCT__
19384a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate"
19399b94acceSBarry Smith /*@C
19409b94acceSBarry Smith    SNESGetSolutionUpdate - Returns the vector where the solution update is
19419b94acceSBarry Smith    stored.
19429b94acceSBarry Smith 
1943c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1944c7afd0dbSLois Curfman McInnes 
19459b94acceSBarry Smith    Input Parameter:
19469b94acceSBarry Smith .  snes - the SNES context
19479b94acceSBarry Smith 
19489b94acceSBarry Smith    Output Parameter:
19499b94acceSBarry Smith .  x - the solution update
19509b94acceSBarry Smith 
195136851e7fSLois Curfman McInnes    Level: advanced
195236851e7fSLois Curfman McInnes 
19539b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update
19549b94acceSBarry Smith 
19559b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction
19569b94acceSBarry Smith @*/
1957dfbe8321SBarry Smith PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x)
19589b94acceSBarry Smith {
19593a40ed3dSBarry Smith   PetscFunctionBegin;
19604482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19614482741eSBarry Smith   PetscValidPointer(x,2);
19629b94acceSBarry Smith   *x = snes->vec_sol_update_always;
19633a40ed3dSBarry Smith   PetscFunctionReturn(0);
19649b94acceSBarry Smith }
19659b94acceSBarry Smith 
19664a2ae208SSatish Balay #undef __FUNCT__
19674a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction"
19689b94acceSBarry Smith /*@C
19693638b69dSLois Curfman McInnes    SNESGetFunction - Returns the vector where the function is stored.
19709b94acceSBarry Smith 
1971c7afd0dbSLois Curfman McInnes    Not Collective, but Vec is parallel if SNES is parallel
1972c7afd0dbSLois Curfman McInnes 
19739b94acceSBarry Smith    Input Parameter:
19749b94acceSBarry Smith .  snes - the SNES context
19759b94acceSBarry Smith 
19769b94acceSBarry Smith    Output Parameter:
19777bf4e008SBarry Smith +  r - the function (or PETSC_NULL)
197800036973SBarry Smith .  ctx - the function context (or PETSC_NULL)
197900036973SBarry Smith -  func - the function (or PETSC_NULL)
19809b94acceSBarry Smith 
198136851e7fSLois Curfman McInnes    Level: advanced
198236851e7fSLois Curfman McInnes 
1983a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function
19849b94acceSBarry Smith 
19854b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution()
19869b94acceSBarry Smith @*/
1987*6849ba73SBarry Smith PetscErrorCode SNESGetFunction(SNES snes,Vec *r,void **ctx,PetscErrorCode (**func)(SNES,Vec,Vec,void*))
19889b94acceSBarry Smith {
19893a40ed3dSBarry Smith   PetscFunctionBegin;
19904482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
19917bf4e008SBarry Smith   if (r)    *r    = snes->vec_func_always;
19927bf4e008SBarry Smith   if (ctx)  *ctx  = snes->funP;
199300036973SBarry Smith   if (func) *func = snes->computefunction;
19943a40ed3dSBarry Smith   PetscFunctionReturn(0);
19959b94acceSBarry Smith }
19969b94acceSBarry Smith 
19974a2ae208SSatish Balay #undef __FUNCT__
19984a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix"
19993c7409f5SSatish Balay /*@C
20003c7409f5SSatish Balay    SNESSetOptionsPrefix - Sets the prefix used for searching for all
2001d850072dSLois Curfman McInnes    SNES options in the database.
20023c7409f5SSatish Balay 
2003fee21e36SBarry Smith    Collective on SNES
2004fee21e36SBarry Smith 
2005c7afd0dbSLois Curfman McInnes    Input Parameter:
2006c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2007c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2008c7afd0dbSLois Curfman McInnes 
2009d850072dSLois Curfman McInnes    Notes:
2010a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2011c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2012d850072dSLois Curfman McInnes 
201336851e7fSLois Curfman McInnes    Level: advanced
201436851e7fSLois Curfman McInnes 
20153c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database
2016a86d99e1SLois Curfman McInnes 
2017a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions()
20183c7409f5SSatish Balay @*/
2019dfbe8321SBarry Smith PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[])
20203c7409f5SSatish Balay {
2021dfbe8321SBarry Smith   PetscErrorCode ierr;
20223c7409f5SSatish Balay 
20233a40ed3dSBarry Smith   PetscFunctionBegin;
20244482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2025639f9d9dSBarry Smith   ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
202694b7f48cSBarry Smith   ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20273a40ed3dSBarry Smith   PetscFunctionReturn(0);
20283c7409f5SSatish Balay }
20293c7409f5SSatish Balay 
20304a2ae208SSatish Balay #undef __FUNCT__
20314a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix"
20323c7409f5SSatish Balay /*@C
2033f525115eSLois Curfman McInnes    SNESAppendOptionsPrefix - Appends to the prefix used for searching for all
2034d850072dSLois Curfman McInnes    SNES options in the database.
20353c7409f5SSatish Balay 
2036fee21e36SBarry Smith    Collective on SNES
2037fee21e36SBarry Smith 
2038c7afd0dbSLois Curfman McInnes    Input Parameters:
2039c7afd0dbSLois Curfman McInnes +  snes - the SNES context
2040c7afd0dbSLois Curfman McInnes -  prefix - the prefix to prepend to all option names
2041c7afd0dbSLois Curfman McInnes 
2042d850072dSLois Curfman McInnes    Notes:
2043a83b1b31SSatish Balay    A hyphen (-) must NOT be given at the beginning of the prefix name.
2044c7afd0dbSLois Curfman McInnes    The first character of all runtime options is AUTOMATICALLY the hyphen.
2045d850072dSLois Curfman McInnes 
204636851e7fSLois Curfman McInnes    Level: advanced
204736851e7fSLois Curfman McInnes 
20483c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database
2049a86d99e1SLois Curfman McInnes 
2050a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix()
20513c7409f5SSatish Balay @*/
2052dfbe8321SBarry Smith PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[])
20533c7409f5SSatish Balay {
2054dfbe8321SBarry Smith   PetscErrorCode ierr;
20553c7409f5SSatish Balay 
20563a40ed3dSBarry Smith   PetscFunctionBegin;
20574482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2058639f9d9dSBarry Smith   ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
205994b7f48cSBarry Smith   ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr);
20603a40ed3dSBarry Smith   PetscFunctionReturn(0);
20613c7409f5SSatish Balay }
20623c7409f5SSatish Balay 
20634a2ae208SSatish Balay #undef __FUNCT__
20644a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix"
20659ab63eb5SSatish Balay /*@C
20663c7409f5SSatish Balay    SNESGetOptionsPrefix - Sets the prefix used for searching for all
20673c7409f5SSatish Balay    SNES options in the database.
20683c7409f5SSatish Balay 
2069c7afd0dbSLois Curfman McInnes    Not Collective
2070c7afd0dbSLois Curfman McInnes 
20713c7409f5SSatish Balay    Input Parameter:
20723c7409f5SSatish Balay .  snes - the SNES context
20733c7409f5SSatish Balay 
20743c7409f5SSatish Balay    Output Parameter:
20753c7409f5SSatish Balay .  prefix - pointer to the prefix string used
20763c7409f5SSatish Balay 
20779ab63eb5SSatish Balay    Notes: On the fortran side, the user should pass in a string 'prifix' of
20789ab63eb5SSatish Balay    sufficient length to hold the prefix.
20799ab63eb5SSatish Balay 
208036851e7fSLois Curfman McInnes    Level: advanced
208136851e7fSLois Curfman McInnes 
20823c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database
2083a86d99e1SLois Curfman McInnes 
2084a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix()
20853c7409f5SSatish Balay @*/
2086dfbe8321SBarry Smith PetscErrorCode SNESGetOptionsPrefix(SNES snes,char *prefix[])
20873c7409f5SSatish Balay {
2088dfbe8321SBarry Smith   PetscErrorCode ierr;
20893c7409f5SSatish Balay 
20903a40ed3dSBarry Smith   PetscFunctionBegin;
20914482741eSBarry Smith   PetscValidHeaderSpecific(snes,SNES_COOKIE,1);
2092639f9d9dSBarry Smith   ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr);
20933a40ed3dSBarry Smith   PetscFunctionReturn(0);
20943c7409f5SSatish Balay }
20953c7409f5SSatish Balay 
2096b2002411SLois Curfman McInnes 
20974a2ae208SSatish Balay #undef __FUNCT__
20984a2ae208SSatish Balay #define __FUNCT__ "SNESRegister"
20993cea93caSBarry Smith /*@C
21003cea93caSBarry Smith   SNESRegister - See SNESRegisterDynamic()
21013cea93caSBarry Smith 
21027f6c08e0SMatthew Knepley   Level: advanced
21033cea93caSBarry Smith @*/
2104*6849ba73SBarry Smith PetscErrorCode SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES))
2105b2002411SLois Curfman McInnes {
2106e2d1d2b7SBarry Smith   char fullname[PETSC_MAX_PATH_LEN];
2107dfbe8321SBarry Smith   PetscErrorCode ierr;
2108b2002411SLois Curfman McInnes 
2109b2002411SLois Curfman McInnes   PetscFunctionBegin;
2110b0a32e0cSBarry Smith   ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr);
2111c134de8dSSatish Balay   ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr);
2112b2002411SLois Curfman McInnes   PetscFunctionReturn(0);
2113b2002411SLois Curfman McInnes }
2114da9b6338SBarry Smith 
2115da9b6338SBarry Smith #undef __FUNCT__
2116da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin"
2117dfbe8321SBarry Smith PetscErrorCode SNESTestLocalMin(SNES snes)
2118da9b6338SBarry Smith {
2119dfbe8321SBarry Smith   PetscErrorCode ierr;
2120dfbe8321SBarry Smith   int N,i,j;
2121da9b6338SBarry Smith   Vec         u,uh,fh;
2122da9b6338SBarry Smith   PetscScalar value;
2123da9b6338SBarry Smith   PetscReal   norm;
2124da9b6338SBarry Smith 
2125da9b6338SBarry Smith   PetscFunctionBegin;
2126da9b6338SBarry Smith   ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr);
2127da9b6338SBarry Smith   ierr = VecDuplicate(u,&uh);CHKERRQ(ierr);
2128da9b6338SBarry Smith   ierr = VecDuplicate(u,&fh);CHKERRQ(ierr);
2129da9b6338SBarry Smith 
2130da9b6338SBarry Smith   /* currently only works for sequential */
2131da9b6338SBarry Smith   ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n");
2132da9b6338SBarry Smith   ierr = VecGetSize(u,&N);CHKERRQ(ierr);
2133da9b6338SBarry Smith   for (i=0; i<N; i++) {
2134da9b6338SBarry Smith     ierr = VecCopy(u,uh);CHKERRQ(ierr);
2135da9b6338SBarry Smith     ierr  = PetscPrintf(PETSC_COMM_WORLD,"i = %d\n",i);CHKERRQ(ierr);
2136da9b6338SBarry Smith     for (j=-10; j<11; j++) {
2137ccae9161SBarry Smith       value = PetscSign(j)*exp(PetscAbs(j)-10.0);
2138da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
21393ab0aad5SBarry Smith       ierr  = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr);
2140da9b6338SBarry Smith       ierr  = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr);
2141da9b6338SBarry Smith       ierr  = PetscPrintf(PETSC_COMM_WORLD,"       j norm %d %18.16e\n",j,norm);CHKERRQ(ierr);
2142da9b6338SBarry Smith       value = -value;
2143da9b6338SBarry Smith       ierr  = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr);
2144da9b6338SBarry Smith     }
2145da9b6338SBarry Smith   }
2146da9b6338SBarry Smith   ierr = VecDestroy(uh);CHKERRQ(ierr);
2147da9b6338SBarry Smith   ierr = VecDestroy(fh);CHKERRQ(ierr);
2148da9b6338SBarry Smith   PetscFunctionReturn(0);
2149da9b6338SBarry Smith }
2150