163dd3a1aSKris Buschelman #define PETSCSNES_DLL 29b94acceSBarry Smith 3e090d566SSatish Balay #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 49b94acceSBarry Smith 54c49b128SBarry Smith PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 68ba1e511SMatthew Knepley PetscFList SNESList = PETSC_NULL; 78ba1e511SMatthew Knepley 88ba1e511SMatthew Knepley /* Logging support */ 963dd3a1aSKris Buschelman PetscCookie PETSCSNES_DLLEXPORT SNES_COOKIE = 0; 106849ba73SBarry Smith PetscEvent SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0; 11a09944afSBarry Smith 12a09944afSBarry Smith #undef __FUNCT__ 134a2ae208SSatish Balay #define __FUNCT__ "SNESView" 147e2c5f70SBarry Smith /*@C 159b94acceSBarry Smith SNESView - Prints the SNES data structure. 169b94acceSBarry Smith 174c49b128SBarry Smith Collective on SNES 18fee21e36SBarry Smith 19c7afd0dbSLois Curfman McInnes Input Parameters: 20c7afd0dbSLois Curfman McInnes + SNES - the SNES context 21c7afd0dbSLois Curfman McInnes - viewer - visualization context 22c7afd0dbSLois Curfman McInnes 239b94acceSBarry Smith Options Database Key: 24c8a8ba5cSLois Curfman McInnes . -snes_view - Calls SNESView() at end of SNESSolve() 259b94acceSBarry Smith 269b94acceSBarry Smith Notes: 279b94acceSBarry Smith The available visualization contexts include 28b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 29b0a32e0cSBarry Smith - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 30c8a8ba5cSLois Curfman McInnes output where only the first processor opens 31c8a8ba5cSLois Curfman McInnes the file. All other processors send their 32c8a8ba5cSLois Curfman McInnes data to the first processor to print. 339b94acceSBarry Smith 343e081fefSLois Curfman McInnes The user can open an alternative visualization context with 35b0a32e0cSBarry Smith PetscViewerASCIIOpen() - output to a specified file. 369b94acceSBarry Smith 3736851e7fSLois Curfman McInnes Level: beginner 3836851e7fSLois Curfman McInnes 399b94acceSBarry Smith .keywords: SNES, view 409b94acceSBarry Smith 41b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen() 429b94acceSBarry Smith @*/ 4363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer) 449b94acceSBarry Smith { 459b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 46dfbe8321SBarry Smith PetscErrorCode ierr; 4794b7f48cSBarry Smith KSP ksp; 48e060cb09SBarry Smith SNESType type; 4932077d6dSBarry Smith PetscTruth iascii,isstring; 509b94acceSBarry Smith 513a40ed3dSBarry Smith PetscFunctionBegin; 524482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 53b0a32e0cSBarry Smith if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm); 544482741eSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 55c9780b6fSBarry Smith PetscCheckSameComm(snes,1,viewer,2); 5674679c65SBarry Smith 5732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 58b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 5932077d6dSBarry Smith if (iascii) { 603a7fca6bSBarry Smith if (snes->prefix) { 613a7fca6bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr); 623a7fca6bSBarry Smith } else { 63b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 643a7fca6bSBarry Smith } 65454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 66454a90a3SBarry Smith if (type) { 67b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 68184914b5SBarry Smith } else { 69b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 70184914b5SBarry Smith } 710ef38995SBarry Smith if (snes->view) { 72b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 730ef38995SBarry Smith ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 74b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 750ef38995SBarry Smith } 7677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 7777d8c4bbSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n", 7870441072SBarry Smith snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr); 7977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 8077431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 819b94acceSBarry Smith if (snes->ksp_ewconv) { 829b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 839b94acceSBarry Smith if (kctx) { 8477431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 85b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 86b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 879b94acceSBarry Smith } 889b94acceSBarry Smith } 890f5bd95cSBarry Smith } else if (isstring) { 90454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 91b0a32e0cSBarry Smith ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 9219bcc07fSBarry Smith } 9394b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 94b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 9594b7f48cSBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 96b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 973a40ed3dSBarry Smith PetscFunctionReturn(0); 989b94acceSBarry Smith } 999b94acceSBarry Smith 10076b2cf59SMatthew Knepley /* 10176b2cf59SMatthew Knepley We retain a list of functions that also take SNES command 10276b2cf59SMatthew Knepley line options. These are called at the end SNESSetFromOptions() 10376b2cf59SMatthew Knepley */ 10476b2cf59SMatthew Knepley #define MAXSETFROMOPTIONS 5 105a7cc72afSBarry Smith static PetscInt numberofsetfromoptions; 1066849ba73SBarry Smith static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 10776b2cf59SMatthew Knepley 108e74ef692SMatthew Knepley #undef __FUNCT__ 109e74ef692SMatthew Knepley #define __FUNCT__ "SNESAddOptionsChecker" 110ac226902SBarry Smith /*@C 11176b2cf59SMatthew Knepley SNESAddOptionsChecker - Adds an additional function to check for SNES options. 11276b2cf59SMatthew Knepley 11376b2cf59SMatthew Knepley Not Collective 11476b2cf59SMatthew Knepley 11576b2cf59SMatthew Knepley Input Parameter: 11676b2cf59SMatthew Knepley . snescheck - function that checks for options 11776b2cf59SMatthew Knepley 11876b2cf59SMatthew Knepley Level: developer 11976b2cf59SMatthew Knepley 12076b2cf59SMatthew Knepley .seealso: SNESSetFromOptions() 12176b2cf59SMatthew Knepley @*/ 12263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 12376b2cf59SMatthew Knepley { 12476b2cf59SMatthew Knepley PetscFunctionBegin; 12576b2cf59SMatthew Knepley if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 12677431f27SBarry Smith SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 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 14670441072SBarry Smith . -snes_atol <abstol> - 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 1522492ecdbSBarry Smith . -snes_no_convergence_test - skip convergence test in nonlinear 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 @*/ 18463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT 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; 19077431f27SBarry Smith PetscInt 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); 21470441072SBarry Smith ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,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 2257aaed0d8SBarry Smith ierr = PetscOptionsTruth("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);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); 25763ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"));CHKERRQ(ierr); 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 @*/ 29563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT 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 @*/ 32363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT 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 @*/ 36363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* 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 @*/ 39263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT 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 @*/ 42263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberUnsuccessfulSteps(SNES snes,PetscInt* 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 @*/ 44763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaximumUnsuccessfulSteps(SNES snes, PetscInt 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 @*/ 47363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaximumUnsuccessfulSteps(SNES snes, PetscInt *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 @*/ 50363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberLinearIterations(SNES snes,PetscInt* 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 @*/ 53663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT 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" 5476849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 548e24b481bSBarry Smith { 549e24b481bSBarry Smith PetscFunctionBegin; 550e24b481bSBarry Smith PetscFunctionReturn(0); 551e24b481bSBarry Smith } 552e24b481bSBarry Smith 5539b94acceSBarry Smith /* -----------------------------------------------------------*/ 5544a2ae208SSatish Balay #undef __FUNCT__ 5554a2ae208SSatish Balay #define __FUNCT__ "SNESCreate" 5569b94acceSBarry Smith /*@C 5579b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5589b94acceSBarry Smith 559c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 560c7afd0dbSLois Curfman McInnes 561c7afd0dbSLois Curfman McInnes Input Parameters: 562c7afd0dbSLois Curfman McInnes + comm - MPI communicator 5639b94acceSBarry Smith 5649b94acceSBarry Smith Output Parameter: 5659b94acceSBarry Smith . outsnes - the new SNES context 5669b94acceSBarry Smith 567c7afd0dbSLois Curfman McInnes Options Database Keys: 568c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 569c7afd0dbSLois Curfman McInnes and no preconditioning matrix 570c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 571c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 572c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 573c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 574c1f60f51SBarry Smith 57536851e7fSLois Curfman McInnes Level: beginner 57636851e7fSLois Curfman McInnes 5779b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5789b94acceSBarry Smith 5794b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES 5809b94acceSBarry Smith @*/ 58163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes) 5829b94acceSBarry Smith { 583dfbe8321SBarry Smith PetscErrorCode ierr; 5849b94acceSBarry Smith SNES snes; 5859b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 58637fcc0dbSBarry Smith 5873a40ed3dSBarry Smith PetscFunctionBegin; 588ed1caa07SMatthew Knepley PetscValidPointer(outsnes,2); 5898ba1e511SMatthew Knepley *outsnes = PETSC_NULL; 5908ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 5918ba1e511SMatthew Knepley ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 5928ba1e511SMatthew Knepley #endif 5938ba1e511SMatthew Knepley 59452e6d16bSBarry Smith ierr = PetscHeaderCreate(snes,_p_SNES,PetscInt,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 595e24b481bSBarry Smith snes->bops->publish = SNESPublish_Petsc; 5969b94acceSBarry Smith snes->max_its = 50; 5979750a799SBarry Smith snes->max_funcs = 10000; 5989b94acceSBarry Smith snes->norm = 0.0; 599b4874afaSBarry Smith snes->rtol = 1.e-8; 600b4874afaSBarry Smith snes->ttol = 0.0; 60170441072SBarry Smith snes->abstol = 1.e-50; 6029b94acceSBarry Smith snes->xtol = 1.e-8; 6034b27c08aSLois Curfman McInnes snes->deltatol = 1.e-12; 6049b94acceSBarry Smith snes->nfuncs = 0; 60550ffb88aSMatthew Knepley snes->numFailures = 0; 60650ffb88aSMatthew Knepley snes->maxFailures = 1; 6077a00f4a9SLois Curfman McInnes snes->linear_its = 0; 608639f9d9dSBarry Smith snes->numbermonitors = 0; 6099b94acceSBarry Smith snes->data = 0; 6109b94acceSBarry Smith snes->view = 0; 6114dc4c822SBarry Smith snes->setupcalled = PETSC_FALSE; 612186905e3SBarry Smith snes->ksp_ewconv = PETSC_FALSE; 6136f24a144SLois Curfman McInnes snes->vwork = 0; 6146f24a144SLois Curfman McInnes snes->nwork = 0; 615758f92a0SBarry Smith snes->conv_hist_len = 0; 616758f92a0SBarry Smith snes->conv_hist_max = 0; 617758f92a0SBarry Smith snes->conv_hist = PETSC_NULL; 618758f92a0SBarry Smith snes->conv_hist_its = PETSC_NULL; 619758f92a0SBarry Smith snes->conv_hist_reset = PETSC_TRUE; 620184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 6219b94acceSBarry Smith 6229b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 623b0a32e0cSBarry Smith ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 62452e6d16bSBarry Smith ierr = PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));CHKERRQ(ierr); 6259b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6269b94acceSBarry Smith kctx->version = 2; 6279b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6289b94acceSBarry Smith this was too large for some test cases */ 6299b94acceSBarry Smith kctx->rtol_last = 0; 6309b94acceSBarry Smith kctx->rtol_max = .9; 6319b94acceSBarry Smith kctx->gamma = 1.0; 6329b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6339b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6349b94acceSBarry Smith kctx->threshold = .1; 6359b94acceSBarry Smith kctx->lresid_last = 0; 6369b94acceSBarry Smith kctx->norm_last = 0; 6379b94acceSBarry Smith 63894b7f48cSBarry Smith ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr); 63952e6d16bSBarry Smith ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr); 6405334005bSBarry Smith 6419b94acceSBarry Smith *outsnes = snes; 64200036973SBarry Smith ierr = PetscPublishAll(snes);CHKERRQ(ierr); 6433a40ed3dSBarry Smith PetscFunctionReturn(0); 6449b94acceSBarry Smith } 6459b94acceSBarry Smith 6464a2ae208SSatish Balay #undef __FUNCT__ 6474a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction" 6489b94acceSBarry Smith /*@C 6499b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6509b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6519b94acceSBarry Smith equations. 6529b94acceSBarry Smith 653fee21e36SBarry Smith Collective on SNES 654fee21e36SBarry Smith 655c7afd0dbSLois Curfman McInnes Input Parameters: 656c7afd0dbSLois Curfman McInnes + snes - the SNES context 657c7afd0dbSLois Curfman McInnes . func - function evaluation routine 658c7afd0dbSLois Curfman McInnes . r - vector to store function value 659c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 660c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6619b94acceSBarry Smith 662c7afd0dbSLois Curfman McInnes Calling sequence of func: 6638d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 664c7afd0dbSLois Curfman McInnes 665313e4042SLois Curfman McInnes . f - function vector 666c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 6679b94acceSBarry Smith 6689b94acceSBarry Smith Notes: 6699b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6709b94acceSBarry Smith $ f'(x) x = -f(x), 671c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 6729b94acceSBarry Smith 67336851e7fSLois Curfman McInnes Level: beginner 67436851e7fSLois Curfman McInnes 6759b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6769b94acceSBarry Smith 677a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6789b94acceSBarry Smith @*/ 67963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 6809b94acceSBarry Smith { 6813a40ed3dSBarry Smith PetscFunctionBegin; 6824482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 6834482741eSBarry Smith PetscValidHeaderSpecific(r,VEC_COOKIE,2); 684c9780b6fSBarry Smith PetscCheckSameComm(snes,1,r,2); 685184914b5SBarry Smith 6869b94acceSBarry Smith snes->computefunction = func; 6879b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6889b94acceSBarry Smith snes->funP = ctx; 6893a40ed3dSBarry Smith PetscFunctionReturn(0); 6909b94acceSBarry Smith } 6919b94acceSBarry Smith 6923ab0aad5SBarry Smith /* --------------------------------------------------------------- */ 6933ab0aad5SBarry Smith #undef __FUNCT__ 6943ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs" 6953ab0aad5SBarry Smith /*@C 6963ab0aad5SBarry Smith SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set 6973ab0aad5SBarry Smith it assumes a zero right hand side. 6983ab0aad5SBarry Smith 6993ab0aad5SBarry Smith Collective on SNES 7003ab0aad5SBarry Smith 7013ab0aad5SBarry Smith Input Parameters: 7023ab0aad5SBarry Smith + snes - the SNES context 7033ab0aad5SBarry Smith - rhs - the right hand side vector or PETSC_NULL for a zero right hand side 7043ab0aad5SBarry Smith 7053ab0aad5SBarry Smith Level: intermediate 7063ab0aad5SBarry Smith 7073ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side 7083ab0aad5SBarry Smith 7091096aae1SMatthew Knepley .seealso: SNESGetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 7103ab0aad5SBarry Smith @*/ 71163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetRhs(SNES snes,Vec rhs) 7123ab0aad5SBarry Smith { 713dfbe8321SBarry Smith PetscErrorCode ierr; 7143ab0aad5SBarry Smith 7153ab0aad5SBarry Smith PetscFunctionBegin; 7163ab0aad5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7173ab0aad5SBarry Smith if (rhs) { 7183ab0aad5SBarry Smith PetscValidHeaderSpecific(rhs,VEC_COOKIE,2); 7193ab0aad5SBarry Smith PetscCheckSameComm(snes,1,rhs,2); 7203ab0aad5SBarry Smith ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr); 7213ab0aad5SBarry Smith } 7223ab0aad5SBarry Smith if (snes->afine) { 7233ab0aad5SBarry Smith ierr = VecDestroy(snes->afine);CHKERRQ(ierr); 7243ab0aad5SBarry Smith } 7253ab0aad5SBarry Smith snes->afine = rhs; 7263ab0aad5SBarry Smith PetscFunctionReturn(0); 7273ab0aad5SBarry Smith } 7283ab0aad5SBarry Smith 7294a2ae208SSatish Balay #undef __FUNCT__ 7301096aae1SMatthew Knepley #define __FUNCT__ "SNESGetRhs" 7311096aae1SMatthew Knepley /*@C 7321096aae1SMatthew Knepley SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 7331096aae1SMatthew Knepley it assumes a zero right hand side. 7341096aae1SMatthew Knepley 7351096aae1SMatthew Knepley Collective on SNES 7361096aae1SMatthew Knepley 7371096aae1SMatthew Knepley Input Parameter: 7381096aae1SMatthew Knepley . snes - the SNES context 7391096aae1SMatthew Knepley 7401096aae1SMatthew Knepley Output Parameter: 7411096aae1SMatthew Knepley . rhs - the right hand side vector or PETSC_NULL for a zero right hand side 7421096aae1SMatthew Knepley 7431096aae1SMatthew Knepley Level: intermediate 7441096aae1SMatthew Knepley 7451096aae1SMatthew Knepley .keywords: SNES, nonlinear, get, function, right hand side 7461096aae1SMatthew Knepley 7471096aae1SMatthew Knepley .seealso: SNESSetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 7481096aae1SMatthew Knepley @*/ 7491096aae1SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs) 7501096aae1SMatthew Knepley { 7511096aae1SMatthew Knepley PetscFunctionBegin; 7521096aae1SMatthew Knepley PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7531096aae1SMatthew Knepley PetscValidPointer(rhs,2); 7541096aae1SMatthew Knepley *rhs = snes->afine; 7551096aae1SMatthew Knepley PetscFunctionReturn(0); 7561096aae1SMatthew Knepley } 7571096aae1SMatthew Knepley 7581096aae1SMatthew Knepley #undef __FUNCT__ 7594a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction" 7609b94acceSBarry Smith /*@ 76136851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7629b94acceSBarry Smith SNESSetFunction(). 7639b94acceSBarry Smith 764c7afd0dbSLois Curfman McInnes Collective on SNES 765c7afd0dbSLois Curfman McInnes 7669b94acceSBarry Smith Input Parameters: 767c7afd0dbSLois Curfman McInnes + snes - the SNES context 768c7afd0dbSLois Curfman McInnes - x - input vector 7699b94acceSBarry Smith 7709b94acceSBarry Smith Output Parameter: 7713638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7729b94acceSBarry Smith 7731bffabb2SLois Curfman McInnes Notes: 77436851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 77536851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 77636851e7fSLois Curfman McInnes themselves. 77736851e7fSLois Curfman McInnes 77836851e7fSLois Curfman McInnes Level: developer 77936851e7fSLois Curfman McInnes 7809b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 7819b94acceSBarry Smith 782a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 7839b94acceSBarry Smith @*/ 78463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y) 7859b94acceSBarry Smith { 786dfbe8321SBarry Smith PetscErrorCode ierr; 7879b94acceSBarry Smith 7883a40ed3dSBarry Smith PetscFunctionBegin; 7894482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7904482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 7914482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,3); 792c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 793c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,3); 794184914b5SBarry Smith 795d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 7961096aae1SMatthew Knepley if (snes->computefunction) { 797d64ed03dSBarry Smith PetscStackPush("SNES user function"); 79819717074SBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); 799d64ed03dSBarry Smith PetscStackPop; 800d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 80119717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 80219717074SBarry Smith } 803d5e45103SBarry Smith CHKERRQ(ierr); 8041096aae1SMatthew Knepley } else if (snes->afine) { 8051096aae1SMatthew Knepley ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 8061096aae1SMatthew Knepley } else { 8071096aae1SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve()."); 8081096aae1SMatthew Knepley } 8093ab0aad5SBarry Smith if (snes->afine) { 8103ab0aad5SBarry Smith PetscScalar mone = -1.0; 8112dcb1b2aSMatthew Knepley ierr = VecAXPY(y,mone,snes->afine);CHKERRQ(ierr); 8123ab0aad5SBarry Smith } 813ae3c334cSLois Curfman McInnes snes->nfuncs++; 814d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 8153a40ed3dSBarry Smith PetscFunctionReturn(0); 8169b94acceSBarry Smith } 8179b94acceSBarry Smith 8184a2ae208SSatish Balay #undef __FUNCT__ 8194a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian" 82062fef451SLois Curfman McInnes /*@ 82162fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 82262fef451SLois Curfman McInnes set with SNESSetJacobian(). 82362fef451SLois Curfman McInnes 824c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 825c7afd0dbSLois Curfman McInnes 82662fef451SLois Curfman McInnes Input Parameters: 827c7afd0dbSLois Curfman McInnes + snes - the SNES context 828c7afd0dbSLois Curfman McInnes - x - input vector 82962fef451SLois Curfman McInnes 83062fef451SLois Curfman McInnes Output Parameters: 831c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 83262fef451SLois Curfman McInnes . B - optional preconditioning matrix 833c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 834fee21e36SBarry Smith 83562fef451SLois Curfman McInnes Notes: 83662fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 83762fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 83862fef451SLois Curfman McInnes 83994b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 840dc5a77f8SLois Curfman McInnes flag parameter. 84162fef451SLois Curfman McInnes 84236851e7fSLois Curfman McInnes Level: developer 84336851e7fSLois Curfman McInnes 84462fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 84562fef451SLois Curfman McInnes 84694b7f48cSBarry Smith .seealso: SNESSetJacobian(), KSPSetOperators() 84762fef451SLois Curfman McInnes @*/ 84863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8499b94acceSBarry Smith { 850dfbe8321SBarry Smith PetscErrorCode ierr; 8513a40ed3dSBarry Smith 8523a40ed3dSBarry Smith PetscFunctionBegin; 8534482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8544482741eSBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE,2); 8554482741eSBarry Smith PetscValidPointer(flg,5); 856c9780b6fSBarry Smith PetscCheckSameComm(snes,1,X,2); 8573a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 858d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 859c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 860d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 8619b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 862d64ed03dSBarry Smith PetscStackPop; 863d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 8646d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 8654482741eSBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 8664482741eSBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 8673a40ed3dSBarry Smith PetscFunctionReturn(0); 8689b94acceSBarry Smith } 8699b94acceSBarry Smith 8704a2ae208SSatish Balay #undef __FUNCT__ 8714a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian" 8729b94acceSBarry Smith /*@C 8739b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 874044dda88SLois Curfman McInnes location to store the matrix. 8759b94acceSBarry Smith 876c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 877c7afd0dbSLois Curfman McInnes 8789b94acceSBarry Smith Input Parameters: 879c7afd0dbSLois Curfman McInnes + snes - the SNES context 8809b94acceSBarry Smith . A - Jacobian matrix 8819b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 8829b94acceSBarry Smith . func - Jacobian evaluation routine 883c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 8842cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 8859b94acceSBarry Smith 8869b94acceSBarry Smith Calling sequence of func: 8878d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 8889b94acceSBarry Smith 889c7afd0dbSLois Curfman McInnes + x - input vector 8909b94acceSBarry Smith . A - Jacobian matrix 8919b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 892ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 89394b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 894c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 8959b94acceSBarry Smith 8969b94acceSBarry Smith Notes: 89794b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 8982cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 899ac21db08SLois Curfman McInnes 900ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9019b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9029b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9039b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9049b94acceSBarry Smith throughout the global iterations. 9059b94acceSBarry Smith 90636851e7fSLois Curfman McInnes Level: beginner 90736851e7fSLois Curfman McInnes 9089b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9099b94acceSBarry Smith 91013f74950SBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), , MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor() 9119b94acceSBarry Smith @*/ 91263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 9139b94acceSBarry Smith { 914dfbe8321SBarry Smith PetscErrorCode ierr; 9153a7fca6bSBarry Smith 9163a40ed3dSBarry Smith PetscFunctionBegin; 9174482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 9184482741eSBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 9194482741eSBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 920c9780b6fSBarry Smith if (A) PetscCheckSameComm(snes,1,A,2); 921c9780b6fSBarry Smith if (B) PetscCheckSameComm(snes,1,B,2); 9223a7fca6bSBarry Smith if (func) snes->computejacobian = func; 9233a7fca6bSBarry Smith if (ctx) snes->jacP = ctx; 9243a7fca6bSBarry Smith if (A) { 9253a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 9269b94acceSBarry Smith snes->jacobian = A; 9273a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 9283a7fca6bSBarry Smith } 9293a7fca6bSBarry Smith if (B) { 9303a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 9319b94acceSBarry Smith snes->jacobian_pre = B; 9323a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 9333a7fca6bSBarry Smith } 93469a612a9SBarry Smith ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 9353a40ed3dSBarry Smith PetscFunctionReturn(0); 9369b94acceSBarry Smith } 93762fef451SLois Curfman McInnes 9384a2ae208SSatish Balay #undef __FUNCT__ 9394a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian" 940c2aafc4cSSatish Balay /*@C 941b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 942b4fd4287SBarry Smith provided context for evaluating the Jacobian. 943b4fd4287SBarry Smith 944c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 945c7afd0dbSLois Curfman McInnes 946b4fd4287SBarry Smith Input Parameter: 947b4fd4287SBarry Smith . snes - the nonlinear solver context 948b4fd4287SBarry Smith 949b4fd4287SBarry Smith Output Parameters: 950c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 951b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 95270e92668SMatthew Knepley . func - location to put Jacobian function (or PETSC_NULL) 95370e92668SMatthew Knepley - ctx - location to stash Jacobian ctx (or PETSC_NULL) 954fee21e36SBarry Smith 95536851e7fSLois Curfman McInnes Level: advanced 95636851e7fSLois Curfman McInnes 957b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 958b4fd4287SBarry Smith @*/ 95970e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 960b4fd4287SBarry Smith { 9613a40ed3dSBarry Smith PetscFunctionBegin; 9624482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 963b4fd4287SBarry Smith if (A) *A = snes->jacobian; 964b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 96500036973SBarry Smith if (func) *func = snes->computejacobian; 96670e92668SMatthew Knepley if (ctx) *ctx = snes->jacP; 9673a40ed3dSBarry Smith PetscFunctionReturn(0); 968b4fd4287SBarry Smith } 969b4fd4287SBarry Smith 9709b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 97163dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 9729b94acceSBarry Smith 9734a2ae208SSatish Balay #undef __FUNCT__ 9744a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp" 9759b94acceSBarry Smith /*@ 9769b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 977272ac6f2SLois Curfman McInnes of a nonlinear solver. 9789b94acceSBarry Smith 979fee21e36SBarry Smith Collective on SNES 980fee21e36SBarry Smith 981c7afd0dbSLois Curfman McInnes Input Parameters: 98270e92668SMatthew Knepley . snes - the SNES context 983c7afd0dbSLois Curfman McInnes 984272ac6f2SLois Curfman McInnes Notes: 985272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 986272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 987272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 988272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 989272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 990272ac6f2SLois Curfman McInnes 99136851e7fSLois Curfman McInnes Level: advanced 99236851e7fSLois Curfman McInnes 9939b94acceSBarry Smith .keywords: SNES, nonlinear, setup 9949b94acceSBarry Smith 9959b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 9969b94acceSBarry Smith @*/ 99770e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 9989b94acceSBarry Smith { 999dfbe8321SBarry Smith PetscErrorCode ierr; 10004b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 10013a40ed3dSBarry Smith 10023a40ed3dSBarry Smith PetscFunctionBegin; 10034482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 10044dc4c822SBarry Smith if (snes->setupcalled) PetscFunctionReturn(0); 10059b94acceSBarry Smith 1006b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1007c1f60f51SBarry Smith /* 1008c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1009dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1010c1f60f51SBarry Smith */ 1011c1f60f51SBarry Smith if (flg) { 1012c1f60f51SBarry Smith Mat J; 10135a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10145a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 101563ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator routines\n"));CHKERRQ(ierr); 10163a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 10173a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 1018c1f60f51SBarry Smith } 101945fc7adcSBarry Smith 1020c8d321e0SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) 102145fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 102245fc7adcSBarry Smith if (flg) { 102345fc7adcSBarry Smith Mat J; 102445fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 102545fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 102645fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 102745fc7adcSBarry Smith } 102832a4b47aSMatthew Knepley #endif 102945fc7adcSBarry Smith 1030b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1031c1f60f51SBarry Smith /* 1032dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1033c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1034c1f60f51SBarry Smith */ 103531615d04SLois Curfman McInnes if (flg) { 1036272ac6f2SLois Curfman McInnes Mat J; 1037b5d62d44SBarry Smith KSP ksp; 103894b7f48cSBarry Smith PC pc; 1039f3ef73ceSBarry Smith 10405a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10413a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 104263ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"));CHKERRQ(ierr); 10433a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 10443a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 10453a7fca6bSBarry Smith 1046f3ef73ceSBarry Smith /* force no preconditioner */ 104794b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1048b5d62d44SBarry Smith ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1049a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1050a9815358SBarry Smith if (!flg) { 1051f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1052272ac6f2SLois Curfman McInnes } 1053a9815358SBarry Smith } 1054f3ef73ceSBarry Smith 10551096aae1SMatthew Knepley if (!snes->vec_func && !snes->afine) { 10561096aae1SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 10571096aae1SMatthew Knepley } 10581096aae1SMatthew Knepley if (!snes->computefunction && !snes->afine) { 10591096aae1SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 10601096aae1SMatthew Knepley } 106129bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1062a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 106329bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1064a8c6a408SBarry Smith } 1065a703fe33SLois Curfman McInnes 1066a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 10674b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 10686831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 106994b7f48cSBarry Smith KSP ksp; 107094b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1071186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1072a703fe33SLois Curfman McInnes } 10734b27c08aSLois Curfman McInnes 1074a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 10757aaed0d8SBarry Smith snes->setupcalled = PETSC_TRUE; 10763a40ed3dSBarry Smith PetscFunctionReturn(0); 10779b94acceSBarry Smith } 10789b94acceSBarry Smith 10794a2ae208SSatish Balay #undef __FUNCT__ 10804a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 10819b94acceSBarry Smith /*@C 10829b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 10839b94acceSBarry Smith with SNESCreate(). 10849b94acceSBarry Smith 1085c7afd0dbSLois Curfman McInnes Collective on SNES 1086c7afd0dbSLois Curfman McInnes 10879b94acceSBarry Smith Input Parameter: 10889b94acceSBarry Smith . snes - the SNES context 10899b94acceSBarry Smith 109036851e7fSLois Curfman McInnes Level: beginner 109136851e7fSLois Curfman McInnes 10929b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 10939b94acceSBarry Smith 109463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 10959b94acceSBarry Smith @*/ 109663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 10979b94acceSBarry Smith { 10986849ba73SBarry Smith PetscErrorCode ierr; 10993a40ed3dSBarry Smith 11003a40ed3dSBarry Smith PetscFunctionBegin; 11014482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 11023a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1103d4bb536fSBarry Smith 1104be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 11050f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1106be0abb6dSBarry Smith 1107e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1108606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 11093a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 11103a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 11113ab0aad5SBarry Smith if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 111294b7f48cSBarry Smith ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1113522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1114*d952e501SBarry Smith ierr = SNESClearMonitor(snes);CHKERRQ(ierr); 1115a79aaaedSSatish Balay ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 11163a40ed3dSBarry Smith PetscFunctionReturn(0); 11179b94acceSBarry Smith } 11189b94acceSBarry Smith 11199b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 11209b94acceSBarry Smith 11214a2ae208SSatish Balay #undef __FUNCT__ 11224a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 11239b94acceSBarry Smith /*@ 1124d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 11259b94acceSBarry Smith 1126c7afd0dbSLois Curfman McInnes Collective on SNES 1127c7afd0dbSLois Curfman McInnes 11289b94acceSBarry Smith Input Parameters: 1129c7afd0dbSLois Curfman McInnes + snes - the SNES context 113070441072SBarry Smith . abstol - absolute convergence tolerance 113133174efeSLois Curfman McInnes . rtol - relative convergence tolerance 113233174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 113333174efeSLois Curfman McInnes of the change in the solution between steps 113433174efeSLois Curfman McInnes . maxit - maximum number of iterations 1135c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1136fee21e36SBarry Smith 113733174efeSLois Curfman McInnes Options Database Keys: 113870441072SBarry Smith + -snes_atol <abstol> - Sets abstol 1139c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1140c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1141c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1142c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 11439b94acceSBarry Smith 1144d7a720efSLois Curfman McInnes Notes: 11459b94acceSBarry Smith The default maximum number of iterations is 50. 11469b94acceSBarry Smith The default maximum number of function evaluations is 1000. 11479b94acceSBarry Smith 114836851e7fSLois Curfman McInnes Level: intermediate 114936851e7fSLois Curfman McInnes 115033174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 11519b94acceSBarry Smith 11522492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance() 11539b94acceSBarry Smith @*/ 115463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 11559b94acceSBarry Smith { 11563a40ed3dSBarry Smith PetscFunctionBegin; 11574482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 115870441072SBarry Smith if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1159d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1160d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1161d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1162d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11633a40ed3dSBarry Smith PetscFunctionReturn(0); 11649b94acceSBarry Smith } 11659b94acceSBarry Smith 11664a2ae208SSatish Balay #undef __FUNCT__ 11674a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11689b94acceSBarry Smith /*@ 116933174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 117033174efeSLois Curfman McInnes 1171c7afd0dbSLois Curfman McInnes Not Collective 1172c7afd0dbSLois Curfman McInnes 117333174efeSLois Curfman McInnes Input Parameters: 1174c7afd0dbSLois Curfman McInnes + snes - the SNES context 117570441072SBarry Smith . abstol - absolute convergence tolerance 117633174efeSLois Curfman McInnes . rtol - relative convergence tolerance 117733174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 117833174efeSLois Curfman McInnes of the change in the solution between steps 117933174efeSLois Curfman McInnes . maxit - maximum number of iterations 1180c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1181fee21e36SBarry Smith 118233174efeSLois Curfman McInnes Notes: 118333174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 118433174efeSLois Curfman McInnes 118536851e7fSLois Curfman McInnes Level: intermediate 118636851e7fSLois Curfman McInnes 118733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 118833174efeSLois Curfman McInnes 118933174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 119033174efeSLois Curfman McInnes @*/ 119163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 119233174efeSLois Curfman McInnes { 11933a40ed3dSBarry Smith PetscFunctionBegin; 11944482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 119570441072SBarry Smith if (abstol) *abstol = snes->abstol; 119633174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 119733174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 119833174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 119933174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 12003a40ed3dSBarry Smith PetscFunctionReturn(0); 120133174efeSLois Curfman McInnes } 120233174efeSLois Curfman McInnes 12034a2ae208SSatish Balay #undef __FUNCT__ 12044a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 120533174efeSLois Curfman McInnes /*@ 12069b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 12079b94acceSBarry Smith 1208fee21e36SBarry Smith Collective on SNES 1209fee21e36SBarry Smith 1210c7afd0dbSLois Curfman McInnes Input Parameters: 1211c7afd0dbSLois Curfman McInnes + snes - the SNES context 1212c7afd0dbSLois Curfman McInnes - tol - tolerance 1213c7afd0dbSLois Curfman McInnes 12149b94acceSBarry Smith Options Database Key: 1215c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 12169b94acceSBarry Smith 121736851e7fSLois Curfman McInnes Level: intermediate 121836851e7fSLois Curfman McInnes 12199b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 12209b94acceSBarry Smith 12212492ecdbSBarry Smith .seealso: SNESSetTolerances() 12229b94acceSBarry Smith @*/ 122363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 12249b94acceSBarry Smith { 12253a40ed3dSBarry Smith PetscFunctionBegin; 12264482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 12279b94acceSBarry Smith snes->deltatol = tol; 12283a40ed3dSBarry Smith PetscFunctionReturn(0); 12299b94acceSBarry Smith } 12309b94acceSBarry Smith 1231df9fa365SBarry Smith /* 1232df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1233df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1234df9fa365SBarry Smith macros instead of functions 1235df9fa365SBarry Smith */ 12364a2ae208SSatish Balay #undef __FUNCT__ 12374a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 123863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1239ce1608b8SBarry Smith { 1240dfbe8321SBarry Smith PetscErrorCode ierr; 1241ce1608b8SBarry Smith 1242ce1608b8SBarry Smith PetscFunctionBegin; 12434482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1244ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1245ce1608b8SBarry Smith PetscFunctionReturn(0); 1246ce1608b8SBarry Smith } 1247ce1608b8SBarry Smith 12484a2ae208SSatish Balay #undef __FUNCT__ 12494a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 125063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1251df9fa365SBarry Smith { 1252dfbe8321SBarry Smith PetscErrorCode ierr; 1253df9fa365SBarry Smith 1254df9fa365SBarry Smith PetscFunctionBegin; 1255df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1256df9fa365SBarry Smith PetscFunctionReturn(0); 1257df9fa365SBarry Smith } 1258df9fa365SBarry Smith 12594a2ae208SSatish Balay #undef __FUNCT__ 12604a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 126163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw) 1262df9fa365SBarry Smith { 1263dfbe8321SBarry Smith PetscErrorCode ierr; 1264df9fa365SBarry Smith 1265df9fa365SBarry Smith PetscFunctionBegin; 1266df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1267df9fa365SBarry Smith PetscFunctionReturn(0); 1268df9fa365SBarry Smith } 1269df9fa365SBarry Smith 12709b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12719b94acceSBarry Smith 12724a2ae208SSatish Balay #undef __FUNCT__ 12734a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12749b94acceSBarry Smith /*@C 12755cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12769b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12779b94acceSBarry Smith progress. 12789b94acceSBarry Smith 1279fee21e36SBarry Smith Collective on SNES 1280fee21e36SBarry Smith 1281c7afd0dbSLois Curfman McInnes Input Parameters: 1282c7afd0dbSLois Curfman McInnes + snes - the SNES context 1283c7afd0dbSLois Curfman McInnes . func - monitoring routine 1284b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1285b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desitre) 1286b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1287b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 12889b94acceSBarry Smith 1289c7afd0dbSLois Curfman McInnes Calling sequence of func: 1290a7cc72afSBarry Smith $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1291c7afd0dbSLois Curfman McInnes 1292c7afd0dbSLois Curfman McInnes + snes - the SNES context 1293c7afd0dbSLois Curfman McInnes . its - iteration number 1294c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 129540a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 12969b94acceSBarry Smith 12979665c990SLois Curfman McInnes Options Database Keys: 1298c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1299c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1300c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1301c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1302c7afd0dbSLois Curfman McInnes been hardwired into a code by 1303c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1304c7afd0dbSLois Curfman McInnes does not cancel those set via 1305c7afd0dbSLois Curfman McInnes the options database. 13069665c990SLois Curfman McInnes 1307639f9d9dSBarry Smith Notes: 13086bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 13096bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 13106bc08f3fSLois Curfman McInnes order in which they were set. 1311639f9d9dSBarry Smith 131236851e7fSLois Curfman McInnes Level: intermediate 131336851e7fSLois Curfman McInnes 13149b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 13159b94acceSBarry Smith 13165cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 13179b94acceSBarry Smith @*/ 131863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 13199b94acceSBarry Smith { 13203a40ed3dSBarry Smith PetscFunctionBegin; 13214482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1322639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 132329bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1324639f9d9dSBarry Smith } 1325639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1326b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1327639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 13283a40ed3dSBarry Smith PetscFunctionReturn(0); 13299b94acceSBarry Smith } 13309b94acceSBarry Smith 13314a2ae208SSatish Balay #undef __FUNCT__ 13324a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 13335cd90555SBarry Smith /*@C 13345cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 13355cd90555SBarry Smith 1336c7afd0dbSLois Curfman McInnes Collective on SNES 1337c7afd0dbSLois Curfman McInnes 13385cd90555SBarry Smith Input Parameters: 13395cd90555SBarry Smith . snes - the SNES context 13405cd90555SBarry Smith 13411a480d89SAdministrator Options Database Key: 1342c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1343c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1344c7afd0dbSLois Curfman McInnes set via the options database 13455cd90555SBarry Smith 13465cd90555SBarry Smith Notes: 13475cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 13485cd90555SBarry Smith 134936851e7fSLois Curfman McInnes Level: intermediate 135036851e7fSLois Curfman McInnes 13515cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 13525cd90555SBarry Smith 13535cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 13545cd90555SBarry Smith @*/ 135563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes) 13565cd90555SBarry Smith { 1357*d952e501SBarry Smith PetscErrorCode ierr; 1358*d952e501SBarry Smith PetscInt i; 1359*d952e501SBarry Smith 13605cd90555SBarry Smith PetscFunctionBegin; 13614482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1362*d952e501SBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1363*d952e501SBarry Smith if (snes->monitordestroy[i]) { 1364*d952e501SBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1365*d952e501SBarry Smith } 1366*d952e501SBarry Smith } 13675cd90555SBarry Smith snes->numbermonitors = 0; 13685cd90555SBarry Smith PetscFunctionReturn(0); 13695cd90555SBarry Smith } 13705cd90555SBarry Smith 13714a2ae208SSatish Balay #undef __FUNCT__ 13724a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13739b94acceSBarry Smith /*@C 13749b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13759b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13769b94acceSBarry Smith 1377fee21e36SBarry Smith Collective on SNES 1378fee21e36SBarry Smith 1379c7afd0dbSLois Curfman McInnes Input Parameters: 1380c7afd0dbSLois Curfman McInnes + snes - the SNES context 1381c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1382c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1383c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 13849b94acceSBarry Smith 1385c7afd0dbSLois Curfman McInnes Calling sequence of func: 13866849ba73SBarry Smith $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1387c7afd0dbSLois Curfman McInnes 1388c7afd0dbSLois Curfman McInnes + snes - the SNES context 1389c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1390184914b5SBarry Smith . reason - reason for convergence/divergence 1391c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 13924b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 13934b27c08aSLois Curfman McInnes - f - 2-norm of function 13949b94acceSBarry Smith 139536851e7fSLois Curfman McInnes Level: advanced 139636851e7fSLois Curfman McInnes 13979b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13989b94acceSBarry Smith 13994b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 14009b94acceSBarry Smith @*/ 140163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 14029b94acceSBarry Smith { 14033a40ed3dSBarry Smith PetscFunctionBegin; 14044482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14059b94acceSBarry Smith (snes)->converged = func; 14069b94acceSBarry Smith (snes)->cnvP = cctx; 14073a40ed3dSBarry Smith PetscFunctionReturn(0); 14089b94acceSBarry Smith } 14099b94acceSBarry Smith 14104a2ae208SSatish Balay #undef __FUNCT__ 14114a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 1412184914b5SBarry Smith /*@C 1413184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1414184914b5SBarry Smith 1415184914b5SBarry Smith Not Collective 1416184914b5SBarry Smith 1417184914b5SBarry Smith Input Parameter: 1418184914b5SBarry Smith . snes - the SNES context 1419184914b5SBarry Smith 1420184914b5SBarry Smith Output Parameter: 1421e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1422184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1423184914b5SBarry Smith 1424184914b5SBarry Smith Level: intermediate 1425184914b5SBarry Smith 1426184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1427184914b5SBarry Smith 1428184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1429184914b5SBarry Smith 14304b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1431184914b5SBarry Smith @*/ 143263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1433184914b5SBarry Smith { 1434184914b5SBarry Smith PetscFunctionBegin; 14354482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14364482741eSBarry Smith PetscValidPointer(reason,2); 1437184914b5SBarry Smith *reason = snes->reason; 1438184914b5SBarry Smith PetscFunctionReturn(0); 1439184914b5SBarry Smith } 1440184914b5SBarry Smith 14414a2ae208SSatish Balay #undef __FUNCT__ 14424a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1443c9005455SLois Curfman McInnes /*@ 1444c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1445c9005455SLois Curfman McInnes 1446fee21e36SBarry Smith Collective on SNES 1447fee21e36SBarry Smith 1448c7afd0dbSLois Curfman McInnes Input Parameters: 1449c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1450c7afd0dbSLois Curfman McInnes . a - array to hold history 1451cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1452758f92a0SBarry Smith . na - size of a and its 145364731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1454758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1455c7afd0dbSLois Curfman McInnes 1456c9005455SLois Curfman McInnes Notes: 14574b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1458c9005455SLois Curfman McInnes at each step. 1459c9005455SLois Curfman McInnes 1460c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1461c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1462c9005455SLois Curfman McInnes during the section of code that is being timed. 1463c9005455SLois Curfman McInnes 146436851e7fSLois Curfman McInnes Level: intermediate 146536851e7fSLois Curfman McInnes 1466c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1467758f92a0SBarry Smith 146808405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1469758f92a0SBarry Smith 1470c9005455SLois Curfman McInnes @*/ 147163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1472c9005455SLois Curfman McInnes { 14733a40ed3dSBarry Smith PetscFunctionBegin; 14744482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14754482741eSBarry Smith if (na) PetscValidScalarPointer(a,2); 1476c9005455SLois Curfman McInnes snes->conv_hist = a; 1477758f92a0SBarry Smith snes->conv_hist_its = its; 1478758f92a0SBarry Smith snes->conv_hist_max = na; 1479758f92a0SBarry Smith snes->conv_hist_reset = reset; 1480758f92a0SBarry Smith PetscFunctionReturn(0); 1481758f92a0SBarry Smith } 1482758f92a0SBarry Smith 14834a2ae208SSatish Balay #undef __FUNCT__ 14844a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 14850c4c9dddSBarry Smith /*@C 1486758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1487758f92a0SBarry Smith 1488758f92a0SBarry Smith Collective on SNES 1489758f92a0SBarry Smith 1490758f92a0SBarry Smith Input Parameter: 1491758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1492758f92a0SBarry Smith 1493758f92a0SBarry Smith Output Parameters: 1494758f92a0SBarry Smith . a - array to hold history 1495758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1496758f92a0SBarry Smith negative if not converged) for each solve. 1497758f92a0SBarry Smith - na - size of a and its 1498758f92a0SBarry Smith 1499758f92a0SBarry Smith Notes: 1500758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1501758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1502758f92a0SBarry Smith 1503758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1504758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1505758f92a0SBarry Smith during the section of code that is being timed. 1506758f92a0SBarry Smith 1507758f92a0SBarry Smith Level: intermediate 1508758f92a0SBarry Smith 1509758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1510758f92a0SBarry Smith 1511758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1512758f92a0SBarry Smith 1513758f92a0SBarry Smith @*/ 151463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1515758f92a0SBarry Smith { 1516758f92a0SBarry Smith PetscFunctionBegin; 15174482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1518758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1519758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1520758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 15213a40ed3dSBarry Smith PetscFunctionReturn(0); 1522c9005455SLois Curfman McInnes } 1523c9005455SLois Curfman McInnes 1524e74ef692SMatthew Knepley #undef __FUNCT__ 1525e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 1526ac226902SBarry Smith /*@C 152776b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 152876b2cf59SMatthew Knepley at the beginning of every step of the iteration. 152976b2cf59SMatthew Knepley 153076b2cf59SMatthew Knepley Collective on SNES 153176b2cf59SMatthew Knepley 153276b2cf59SMatthew Knepley Input Parameters: 153376b2cf59SMatthew Knepley . snes - The nonlinear solver context 153476b2cf59SMatthew Knepley . func - The function 153576b2cf59SMatthew Knepley 153676b2cf59SMatthew Knepley Calling sequence of func: 1537b5d30489SBarry Smith . func (SNES snes, PetscInt step); 153876b2cf59SMatthew Knepley 153976b2cf59SMatthew Knepley . step - The current step of the iteration 154076b2cf59SMatthew Knepley 154176b2cf59SMatthew Knepley Level: intermediate 154276b2cf59SMatthew Knepley 154376b2cf59SMatthew Knepley .keywords: SNES, update 1544b5d30489SBarry Smith 154576b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 154676b2cf59SMatthew Knepley @*/ 154763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 154876b2cf59SMatthew Knepley { 154976b2cf59SMatthew Knepley PetscFunctionBegin; 15504482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 155176b2cf59SMatthew Knepley snes->update = func; 155276b2cf59SMatthew Knepley PetscFunctionReturn(0); 155376b2cf59SMatthew Knepley } 155476b2cf59SMatthew Knepley 1555e74ef692SMatthew Knepley #undef __FUNCT__ 1556e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 155776b2cf59SMatthew Knepley /*@ 155876b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 155976b2cf59SMatthew Knepley 156076b2cf59SMatthew Knepley Not collective 156176b2cf59SMatthew Knepley 156276b2cf59SMatthew Knepley Input Parameters: 156376b2cf59SMatthew Knepley . snes - The nonlinear solver context 156476b2cf59SMatthew Knepley . step - The current step of the iteration 156576b2cf59SMatthew Knepley 1566205452f4SMatthew Knepley Level: intermediate 1567205452f4SMatthew Knepley 156876b2cf59SMatthew Knepley .keywords: SNES, update 156976b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 157076b2cf59SMatthew Knepley @*/ 157163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 157276b2cf59SMatthew Knepley { 157376b2cf59SMatthew Knepley PetscFunctionBegin; 157476b2cf59SMatthew Knepley PetscFunctionReturn(0); 157576b2cf59SMatthew Knepley } 157676b2cf59SMatthew Knepley 15774a2ae208SSatish Balay #undef __FUNCT__ 15784a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 15799b94acceSBarry Smith /* 15809b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 15819b94acceSBarry Smith positive parameter delta. 15829b94acceSBarry Smith 15839b94acceSBarry Smith Input Parameters: 1584c7afd0dbSLois Curfman McInnes + snes - the SNES context 15859b94acceSBarry Smith . y - approximate solution of linear system 15869b94acceSBarry Smith . fnorm - 2-norm of current function 1587c7afd0dbSLois Curfman McInnes - delta - trust region size 15889b94acceSBarry Smith 15899b94acceSBarry Smith Output Parameters: 1590c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 15919b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 15929b94acceSBarry Smith region, and exceeds zero otherwise. 1593c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 15949b94acceSBarry Smith 15959b94acceSBarry Smith Note: 15964b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 15979b94acceSBarry Smith is set to be the maximum allowable step size. 15989b94acceSBarry Smith 15999b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16009b94acceSBarry Smith */ 1601dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 16029b94acceSBarry Smith { 1603064f8208SBarry Smith PetscReal nrm; 1604ea709b57SSatish Balay PetscScalar cnorm; 1605dfbe8321SBarry Smith PetscErrorCode ierr; 16063a40ed3dSBarry Smith 16073a40ed3dSBarry Smith PetscFunctionBegin; 16084482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16094482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1610c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,2); 1611184914b5SBarry Smith 1612064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1613064f8208SBarry Smith if (nrm > *delta) { 1614064f8208SBarry Smith nrm = *delta/nrm; 1615064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1616064f8208SBarry Smith cnorm = nrm; 16172dcb1b2aSMatthew Knepley ierr = VecScale(y,cnorm);CHKERRQ(ierr); 16189b94acceSBarry Smith *ynorm = *delta; 16199b94acceSBarry Smith } else { 16209b94acceSBarry Smith *gpnorm = 0.0; 1621064f8208SBarry Smith *ynorm = nrm; 16229b94acceSBarry Smith } 16233a40ed3dSBarry Smith PetscFunctionReturn(0); 16249b94acceSBarry Smith } 16259b94acceSBarry Smith 16264a2ae208SSatish Balay #undef __FUNCT__ 16274a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 16289b94acceSBarry Smith /*@ 1629f69a0ea3SMatthew Knepley SNESSolve - Solves a nonlinear system F(x) = b. 1630f69a0ea3SMatthew Knepley Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 16319b94acceSBarry Smith 1632c7afd0dbSLois Curfman McInnes Collective on SNES 1633c7afd0dbSLois Curfman McInnes 1634b2002411SLois Curfman McInnes Input Parameters: 1635c7afd0dbSLois Curfman McInnes + snes - the SNES context 1636f69a0ea3SMatthew Knepley . b - the constant part of the equation, or PETSC_NULL to use zero. 163770e92668SMatthew Knepley - x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution() 16389b94acceSBarry Smith 1639b2002411SLois Curfman McInnes Notes: 16408ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 16418ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 16428ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 16438ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 16448ddd3da0SLois Curfman McInnes 164536851e7fSLois Curfman McInnes Level: beginner 164636851e7fSLois Curfman McInnes 16479b94acceSBarry Smith .keywords: SNES, nonlinear, solve 16489b94acceSBarry Smith 164970e92668SMatthew Knepley .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution() 16509b94acceSBarry Smith @*/ 1651f69a0ea3SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 16529b94acceSBarry Smith { 1653dfbe8321SBarry Smith PetscErrorCode ierr; 1654f1af5d2fSBarry Smith PetscTruth flg; 1655052efed2SBarry Smith 16563a40ed3dSBarry Smith PetscFunctionBegin; 16574482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16581302d50aSBarry Smith if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 165974637425SBarry Smith 1660f69a0ea3SMatthew Knepley if (b) { 1661f69a0ea3SMatthew Knepley ierr = SNESSetRhs(snes, b); CHKERRQ(ierr); 16621096aae1SMatthew Knepley if (!snes->vec_func) { 16631096aae1SMatthew Knepley Vec r; 16641096aae1SMatthew Knepley 16651096aae1SMatthew Knepley ierr = VecDuplicate(b, &r); CHKERRQ(ierr); 16661096aae1SMatthew Knepley ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr); 16671096aae1SMatthew Knepley } 1668f69a0ea3SMatthew Knepley } 166970e92668SMatthew Knepley if (x) { 1670f69a0ea3SMatthew Knepley PetscValidHeaderSpecific(x,VEC_COOKIE,3); 1671f69a0ea3SMatthew Knepley PetscCheckSameComm(snes,1,x,3); 167270e92668SMatthew Knepley } else { 167370e92668SMatthew Knepley ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr); 167470e92668SMatthew Knepley if (!x) { 167570e92668SMatthew Knepley ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr); 167670e92668SMatthew Knepley } 167770e92668SMatthew Knepley } 167870e92668SMatthew Knepley snes->vec_sol = snes->vec_sol_always = x; 167970e92668SMatthew Knepley if (!snes->setupcalled) { 168070e92668SMatthew Knepley ierr = SNESSetUp(snes);CHKERRQ(ierr); 168170e92668SMatthew Knepley } 1682abc0a331SBarry Smith if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1683d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 168450ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1685d5e45103SBarry Smith 1686d5e45103SBarry Smith ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN); 1687d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 168838f152feSBarry Smith /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */ 168919717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr); 1690d5e45103SBarry Smith } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) { 1691d5e45103SBarry Smith /* translate exception into SNES not converged reason */ 1692d5e45103SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 169338f152feSBarry Smith ierr = 0; 169438f152feSBarry Smith } 169538f152feSBarry Smith CHKERRQ(ierr); 1696d5e45103SBarry Smith 1697d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1698b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 16998b179ff8SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1700da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1701da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 17025968eb51SBarry Smith if (snes->printreason) { 17035968eb51SBarry Smith if (snes->reason > 0) { 17049dcbbd2bSBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 17055968eb51SBarry Smith } else { 17069dcbbd2bSBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 17075968eb51SBarry Smith } 17085968eb51SBarry Smith } 17095968eb51SBarry Smith 17103a40ed3dSBarry Smith PetscFunctionReturn(0); 17119b94acceSBarry Smith } 17129b94acceSBarry Smith 17139b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 17149b94acceSBarry Smith 17154a2ae208SSatish Balay #undef __FUNCT__ 17164a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 171782bf6240SBarry Smith /*@C 17184b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 17199b94acceSBarry Smith 1720fee21e36SBarry Smith Collective on SNES 1721fee21e36SBarry Smith 1722c7afd0dbSLois Curfman McInnes Input Parameters: 1723c7afd0dbSLois Curfman McInnes + snes - the SNES context 1724454a90a3SBarry Smith - type - a known method 1725c7afd0dbSLois Curfman McInnes 1726c7afd0dbSLois Curfman McInnes Options Database Key: 1727454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1728c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1729ae12b187SLois Curfman McInnes 17309b94acceSBarry Smith Notes: 1731e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 17324b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1733c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17344b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1735c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17369b94acceSBarry Smith 1737ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1738ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1739ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1740ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1741ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1742ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1743ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1744ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1745ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1746b0a32e0cSBarry Smith appropriate method. 174736851e7fSLois Curfman McInnes 174836851e7fSLois Curfman McInnes Level: intermediate 1749a703fe33SLois Curfman McInnes 1750454a90a3SBarry Smith .keywords: SNES, set, type 1751435da068SBarry Smith 1752435da068SBarry Smith .seealso: SNESType, SNESCreate() 1753435da068SBarry Smith 17549b94acceSBarry Smith @*/ 1755e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type) 17569b94acceSBarry Smith { 1757dfbe8321SBarry Smith PetscErrorCode ierr,(*r)(SNES); 17586831982aSBarry Smith PetscTruth match; 17593a40ed3dSBarry Smith 17603a40ed3dSBarry Smith PetscFunctionBegin; 17614482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17624482741eSBarry Smith PetscValidCharPointer(type,2); 176382bf6240SBarry Smith 17646831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 17650f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 176682bf6240SBarry Smith 176782bf6240SBarry Smith if (snes->setupcalled) { 17684dc4c822SBarry Smith snes->setupcalled = PETSC_FALSE; 1769e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 177082bf6240SBarry Smith snes->data = 0; 177182bf6240SBarry Smith } 17727d1a2b2bSBarry Smith 17739b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 177482bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1775b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1776958c9bccSBarry Smith if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1777606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 177882bf6240SBarry Smith snes->data = 0; 17793a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1780454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 17813a40ed3dSBarry Smith PetscFunctionReturn(0); 17829b94acceSBarry Smith } 17839b94acceSBarry Smith 1784a847f771SSatish Balay 17859b94acceSBarry Smith /* --------------------------------------------------------------------- */ 17864a2ae208SSatish Balay #undef __FUNCT__ 17874a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 17889b94acceSBarry Smith /*@C 17899b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1790f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 17919b94acceSBarry Smith 1792fee21e36SBarry Smith Not Collective 1793fee21e36SBarry Smith 179436851e7fSLois Curfman McInnes Level: advanced 179536851e7fSLois Curfman McInnes 17969b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 17979b94acceSBarry Smith 17989b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 17999b94acceSBarry Smith @*/ 180063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 18019b94acceSBarry Smith { 1802dfbe8321SBarry Smith PetscErrorCode ierr; 180382bf6240SBarry Smith 18043a40ed3dSBarry Smith PetscFunctionBegin; 180582bf6240SBarry Smith if (SNESList) { 1806b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 180782bf6240SBarry Smith SNESList = 0; 18089b94acceSBarry Smith } 18094c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 18103a40ed3dSBarry Smith PetscFunctionReturn(0); 18119b94acceSBarry Smith } 18129b94acceSBarry Smith 18134a2ae208SSatish Balay #undef __FUNCT__ 18144a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 18159b94acceSBarry Smith /*@C 18169a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 18179b94acceSBarry Smith 1818c7afd0dbSLois Curfman McInnes Not Collective 1819c7afd0dbSLois Curfman McInnes 18209b94acceSBarry Smith Input Parameter: 18214b0e389bSBarry Smith . snes - nonlinear solver context 18229b94acceSBarry Smith 18239b94acceSBarry Smith Output Parameter: 18243a7fca6bSBarry Smith . type - SNES method (a character string) 18259b94acceSBarry Smith 182636851e7fSLois Curfman McInnes Level: intermediate 182736851e7fSLois Curfman McInnes 1828454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 18299b94acceSBarry Smith @*/ 183063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type) 18319b94acceSBarry Smith { 18323a40ed3dSBarry Smith PetscFunctionBegin; 18334482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18344482741eSBarry Smith PetscValidPointer(type,2); 1835454a90a3SBarry Smith *type = snes->type_name; 18363a40ed3dSBarry Smith PetscFunctionReturn(0); 18379b94acceSBarry Smith } 18389b94acceSBarry Smith 18394a2ae208SSatish Balay #undef __FUNCT__ 18404a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 18419b94acceSBarry Smith /*@C 18429b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 18439b94acceSBarry Smith stored. 18449b94acceSBarry Smith 1845c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1846c7afd0dbSLois Curfman McInnes 18479b94acceSBarry Smith Input Parameter: 18489b94acceSBarry Smith . snes - the SNES context 18499b94acceSBarry Smith 18509b94acceSBarry Smith Output Parameter: 18519b94acceSBarry Smith . x - the solution 18529b94acceSBarry Smith 185370e92668SMatthew Knepley Level: intermediate 185436851e7fSLois Curfman McInnes 18559b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 18569b94acceSBarry Smith 185770e92668SMatthew Knepley .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 18589b94acceSBarry Smith @*/ 185963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 18609b94acceSBarry Smith { 18613a40ed3dSBarry Smith PetscFunctionBegin; 18624482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18634482741eSBarry Smith PetscValidPointer(x,2); 18649b94acceSBarry Smith *x = snes->vec_sol_always; 18653a40ed3dSBarry Smith PetscFunctionReturn(0); 18669b94acceSBarry Smith } 18679b94acceSBarry Smith 18684a2ae208SSatish Balay #undef __FUNCT__ 186970e92668SMatthew Knepley #define __FUNCT__ "SNESSetSolution" 187070e92668SMatthew Knepley /*@C 187170e92668SMatthew Knepley SNESSetSolution - Sets the vector where the approximate solution is stored. 187270e92668SMatthew Knepley 187370e92668SMatthew Knepley Not Collective, but Vec is parallel if SNES is parallel 187470e92668SMatthew Knepley 187570e92668SMatthew Knepley Input Parameters: 187670e92668SMatthew Knepley + snes - the SNES context 187770e92668SMatthew Knepley - x - the solution 187870e92668SMatthew Knepley 187970e92668SMatthew Knepley Output Parameter: 188070e92668SMatthew Knepley 188170e92668SMatthew Knepley Level: intermediate 188270e92668SMatthew Knepley 188370e92668SMatthew Knepley .keywords: SNES, nonlinear, set, solution 188470e92668SMatthew Knepley 188570e92668SMatthew Knepley .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 188670e92668SMatthew Knepley @*/ 188770e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x) 188870e92668SMatthew Knepley { 188970e92668SMatthew Knepley PetscFunctionBegin; 189070e92668SMatthew Knepley PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 189170e92668SMatthew Knepley PetscValidHeaderSpecific(x,VEC_COOKIE,2); 189270e92668SMatthew Knepley PetscCheckSameComm(snes,1,x,2); 189370e92668SMatthew Knepley snes->vec_sol_always = x; 189470e92668SMatthew Knepley PetscFunctionReturn(0); 189570e92668SMatthew Knepley } 189670e92668SMatthew Knepley 189770e92668SMatthew Knepley #undef __FUNCT__ 18984a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 18999b94acceSBarry Smith /*@C 19009b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 19019b94acceSBarry Smith stored. 19029b94acceSBarry Smith 1903c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1904c7afd0dbSLois Curfman McInnes 19059b94acceSBarry Smith Input Parameter: 19069b94acceSBarry Smith . snes - the SNES context 19079b94acceSBarry Smith 19089b94acceSBarry Smith Output Parameter: 19099b94acceSBarry Smith . x - the solution update 19109b94acceSBarry Smith 191136851e7fSLois Curfman McInnes Level: advanced 191236851e7fSLois Curfman McInnes 19139b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 19149b94acceSBarry Smith 19159b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 19169b94acceSBarry Smith @*/ 191763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 19189b94acceSBarry Smith { 19193a40ed3dSBarry Smith PetscFunctionBegin; 19204482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19214482741eSBarry Smith PetscValidPointer(x,2); 19229b94acceSBarry Smith *x = snes->vec_sol_update_always; 19233a40ed3dSBarry Smith PetscFunctionReturn(0); 19249b94acceSBarry Smith } 19259b94acceSBarry Smith 19264a2ae208SSatish Balay #undef __FUNCT__ 19274a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 19289b94acceSBarry Smith /*@C 19293638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 19309b94acceSBarry Smith 1931c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1932c7afd0dbSLois Curfman McInnes 19339b94acceSBarry Smith Input Parameter: 19349b94acceSBarry Smith . snes - the SNES context 19359b94acceSBarry Smith 19369b94acceSBarry Smith Output Parameter: 19377bf4e008SBarry Smith + r - the function (or PETSC_NULL) 193870e92668SMatthew Knepley . func - the function (or PETSC_NULL) 193970e92668SMatthew Knepley - ctx - the function context (or PETSC_NULL) 19409b94acceSBarry Smith 194136851e7fSLois Curfman McInnes Level: advanced 194236851e7fSLois Curfman McInnes 1943a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 19449b94acceSBarry Smith 19454b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 19469b94acceSBarry Smith @*/ 194770e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 19489b94acceSBarry Smith { 19493a40ed3dSBarry Smith PetscFunctionBegin; 19504482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19517bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 195200036973SBarry Smith if (func) *func = snes->computefunction; 195370e92668SMatthew Knepley if (ctx) *ctx = snes->funP; 19543a40ed3dSBarry Smith PetscFunctionReturn(0); 19559b94acceSBarry Smith } 19569b94acceSBarry Smith 19574a2ae208SSatish Balay #undef __FUNCT__ 19584a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 19593c7409f5SSatish Balay /*@C 19603c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1961d850072dSLois Curfman McInnes SNES options in the database. 19623c7409f5SSatish Balay 1963fee21e36SBarry Smith Collective on SNES 1964fee21e36SBarry Smith 1965c7afd0dbSLois Curfman McInnes Input Parameter: 1966c7afd0dbSLois Curfman McInnes + snes - the SNES context 1967c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1968c7afd0dbSLois Curfman McInnes 1969d850072dSLois Curfman McInnes Notes: 1970a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1971c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1972d850072dSLois Curfman McInnes 197336851e7fSLois Curfman McInnes Level: advanced 197436851e7fSLois Curfman McInnes 19753c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1976a86d99e1SLois Curfman McInnes 1977a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19783c7409f5SSatish Balay @*/ 197963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 19803c7409f5SSatish Balay { 1981dfbe8321SBarry Smith PetscErrorCode ierr; 19823c7409f5SSatish Balay 19833a40ed3dSBarry Smith PetscFunctionBegin; 19844482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1985639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 198694b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 19873a40ed3dSBarry Smith PetscFunctionReturn(0); 19883c7409f5SSatish Balay } 19893c7409f5SSatish Balay 19904a2ae208SSatish Balay #undef __FUNCT__ 19914a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 19923c7409f5SSatish Balay /*@C 1993f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1994d850072dSLois Curfman McInnes SNES options in the database. 19953c7409f5SSatish Balay 1996fee21e36SBarry Smith Collective on SNES 1997fee21e36SBarry Smith 1998c7afd0dbSLois Curfman McInnes Input Parameters: 1999c7afd0dbSLois Curfman McInnes + snes - the SNES context 2000c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2001c7afd0dbSLois Curfman McInnes 2002d850072dSLois Curfman McInnes Notes: 2003a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2004c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2005d850072dSLois Curfman McInnes 200636851e7fSLois Curfman McInnes Level: advanced 200736851e7fSLois Curfman McInnes 20083c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2009a86d99e1SLois Curfman McInnes 2010a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20113c7409f5SSatish Balay @*/ 201263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 20133c7409f5SSatish Balay { 2014dfbe8321SBarry Smith PetscErrorCode ierr; 20153c7409f5SSatish Balay 20163a40ed3dSBarry Smith PetscFunctionBegin; 20174482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2018639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 201994b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20203a40ed3dSBarry Smith PetscFunctionReturn(0); 20213c7409f5SSatish Balay } 20223c7409f5SSatish Balay 20234a2ae208SSatish Balay #undef __FUNCT__ 20244a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 20259ab63eb5SSatish Balay /*@C 20263c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20273c7409f5SSatish Balay SNES options in the database. 20283c7409f5SSatish Balay 2029c7afd0dbSLois Curfman McInnes Not Collective 2030c7afd0dbSLois Curfman McInnes 20313c7409f5SSatish Balay Input Parameter: 20323c7409f5SSatish Balay . snes - the SNES context 20333c7409f5SSatish Balay 20343c7409f5SSatish Balay Output Parameter: 20353c7409f5SSatish Balay . prefix - pointer to the prefix string used 20363c7409f5SSatish Balay 20379ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 20389ab63eb5SSatish Balay sufficient length to hold the prefix. 20399ab63eb5SSatish Balay 204036851e7fSLois Curfman McInnes Level: advanced 204136851e7fSLois Curfman McInnes 20423c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2043a86d99e1SLois Curfman McInnes 2044a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20453c7409f5SSatish Balay @*/ 2046e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 20473c7409f5SSatish Balay { 2048dfbe8321SBarry Smith PetscErrorCode ierr; 20493c7409f5SSatish Balay 20503a40ed3dSBarry Smith PetscFunctionBegin; 20514482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2052639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 20533a40ed3dSBarry Smith PetscFunctionReturn(0); 20543c7409f5SSatish Balay } 20553c7409f5SSatish Balay 2056b2002411SLois Curfman McInnes 20574a2ae208SSatish Balay #undef __FUNCT__ 20584a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 20593cea93caSBarry Smith /*@C 20603cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 20613cea93caSBarry Smith 20627f6c08e0SMatthew Knepley Level: advanced 20633cea93caSBarry Smith @*/ 206463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2065b2002411SLois Curfman McInnes { 2066e2d1d2b7SBarry Smith char fullname[PETSC_MAX_PATH_LEN]; 2067dfbe8321SBarry Smith PetscErrorCode ierr; 2068b2002411SLois Curfman McInnes 2069b2002411SLois Curfman McInnes PetscFunctionBegin; 2070b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2071c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2072b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2073b2002411SLois Curfman McInnes } 2074da9b6338SBarry Smith 2075da9b6338SBarry Smith #undef __FUNCT__ 2076da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 207763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2078da9b6338SBarry Smith { 2079dfbe8321SBarry Smith PetscErrorCode ierr; 208077431f27SBarry Smith PetscInt N,i,j; 2081da9b6338SBarry Smith Vec u,uh,fh; 2082da9b6338SBarry Smith PetscScalar value; 2083da9b6338SBarry Smith PetscReal norm; 2084da9b6338SBarry Smith 2085da9b6338SBarry Smith PetscFunctionBegin; 2086da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2087da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2088da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2089da9b6338SBarry Smith 2090da9b6338SBarry Smith /* currently only works for sequential */ 2091da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2092da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2093da9b6338SBarry Smith for (i=0; i<N; i++) { 2094da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 209577431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2096da9b6338SBarry Smith for (j=-10; j<11; j++) { 2097ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2098da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 20993ab0aad5SBarry Smith ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2100da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 210177431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2102da9b6338SBarry Smith value = -value; 2103da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2104da9b6338SBarry Smith } 2105da9b6338SBarry Smith } 2106da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2107da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2108da9b6338SBarry Smith PetscFunctionReturn(0); 2109da9b6338SBarry Smith } 2110