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"); 798e9a2bbcdSBarry Smith CHKMEMQ; 79919717074SBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); 800e9a2bbcdSBarry Smith CHKMEMQ; 801d64ed03dSBarry Smith PetscStackPop; 802d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 80319717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 80419717074SBarry Smith } 805d5e45103SBarry Smith CHKERRQ(ierr); 8061096aae1SMatthew Knepley } else if (snes->afine) { 8071096aae1SMatthew Knepley ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 8081096aae1SMatthew Knepley } else { 8091096aae1SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve()."); 8101096aae1SMatthew Knepley } 8113ab0aad5SBarry Smith if (snes->afine) { 812*016dedfbSBarry Smith ierr = VecAXPY(y,-1.0,snes->afine);CHKERRQ(ierr); 8133ab0aad5SBarry Smith } 814ae3c334cSLois Curfman McInnes snes->nfuncs++; 815d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 8163a40ed3dSBarry Smith PetscFunctionReturn(0); 8179b94acceSBarry Smith } 8189b94acceSBarry Smith 8194a2ae208SSatish Balay #undef __FUNCT__ 8204a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian" 82162fef451SLois Curfman McInnes /*@ 82262fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 82362fef451SLois Curfman McInnes set with SNESSetJacobian(). 82462fef451SLois Curfman McInnes 825c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 826c7afd0dbSLois Curfman McInnes 82762fef451SLois Curfman McInnes Input Parameters: 828c7afd0dbSLois Curfman McInnes + snes - the SNES context 829c7afd0dbSLois Curfman McInnes - x - input vector 83062fef451SLois Curfman McInnes 83162fef451SLois Curfman McInnes Output Parameters: 832c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 83362fef451SLois Curfman McInnes . B - optional preconditioning matrix 834c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 835fee21e36SBarry Smith 83662fef451SLois Curfman McInnes Notes: 83762fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 83862fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 83962fef451SLois Curfman McInnes 84094b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 841dc5a77f8SLois Curfman McInnes flag parameter. 84262fef451SLois Curfman McInnes 84336851e7fSLois Curfman McInnes Level: developer 84436851e7fSLois Curfman McInnes 84562fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 84662fef451SLois Curfman McInnes 84794b7f48cSBarry Smith .seealso: SNESSetJacobian(), KSPSetOperators() 84862fef451SLois Curfman McInnes @*/ 84963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8509b94acceSBarry Smith { 851dfbe8321SBarry Smith PetscErrorCode ierr; 8523a40ed3dSBarry Smith 8533a40ed3dSBarry Smith PetscFunctionBegin; 8544482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8554482741eSBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE,2); 8564482741eSBarry Smith PetscValidPointer(flg,5); 857c9780b6fSBarry Smith PetscCheckSameComm(snes,1,X,2); 8583a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 859d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 860c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 861d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 8629b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 863d64ed03dSBarry Smith PetscStackPop; 864d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 8656d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 8664482741eSBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 8674482741eSBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 8683a40ed3dSBarry Smith PetscFunctionReturn(0); 8699b94acceSBarry Smith } 8709b94acceSBarry Smith 8714a2ae208SSatish Balay #undef __FUNCT__ 8724a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian" 8739b94acceSBarry Smith /*@C 8749b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 875044dda88SLois Curfman McInnes location to store the matrix. 8769b94acceSBarry Smith 877c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 878c7afd0dbSLois Curfman McInnes 8799b94acceSBarry Smith Input Parameters: 880c7afd0dbSLois Curfman McInnes + snes - the SNES context 8819b94acceSBarry Smith . A - Jacobian matrix 8829b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 8839b94acceSBarry Smith . func - Jacobian evaluation routine 884c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 8852cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 8869b94acceSBarry Smith 8879b94acceSBarry Smith Calling sequence of func: 8888d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 8899b94acceSBarry Smith 890c7afd0dbSLois Curfman McInnes + x - input vector 8919b94acceSBarry Smith . A - Jacobian matrix 8929b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 893ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 89494b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 895c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 8969b94acceSBarry Smith 8979b94acceSBarry Smith Notes: 89894b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 8992cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 900ac21db08SLois Curfman McInnes 901ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9029b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9039b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9049b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9059b94acceSBarry Smith throughout the global iterations. 9069b94acceSBarry Smith 90736851e7fSLois Curfman McInnes Level: beginner 90836851e7fSLois Curfman McInnes 9099b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9109b94acceSBarry Smith 91113f74950SBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), , MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor() 9129b94acceSBarry Smith @*/ 91363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 9149b94acceSBarry Smith { 915dfbe8321SBarry Smith PetscErrorCode ierr; 9163a7fca6bSBarry Smith 9173a40ed3dSBarry Smith PetscFunctionBegin; 9184482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 9194482741eSBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 9204482741eSBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 921c9780b6fSBarry Smith if (A) PetscCheckSameComm(snes,1,A,2); 922c9780b6fSBarry Smith if (B) PetscCheckSameComm(snes,1,B,2); 9233a7fca6bSBarry Smith if (func) snes->computejacobian = func; 9243a7fca6bSBarry Smith if (ctx) snes->jacP = ctx; 9253a7fca6bSBarry Smith if (A) { 9263a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 9279b94acceSBarry Smith snes->jacobian = A; 9283a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 9293a7fca6bSBarry Smith } 9303a7fca6bSBarry Smith if (B) { 9313a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 9329b94acceSBarry Smith snes->jacobian_pre = B; 9333a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 9343a7fca6bSBarry Smith } 93569a612a9SBarry Smith ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 9363a40ed3dSBarry Smith PetscFunctionReturn(0); 9379b94acceSBarry Smith } 93862fef451SLois Curfman McInnes 9394a2ae208SSatish Balay #undef __FUNCT__ 9404a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian" 941c2aafc4cSSatish Balay /*@C 942b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 943b4fd4287SBarry Smith provided context for evaluating the Jacobian. 944b4fd4287SBarry Smith 945c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 946c7afd0dbSLois Curfman McInnes 947b4fd4287SBarry Smith Input Parameter: 948b4fd4287SBarry Smith . snes - the nonlinear solver context 949b4fd4287SBarry Smith 950b4fd4287SBarry Smith Output Parameters: 951c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 952b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 95370e92668SMatthew Knepley . func - location to put Jacobian function (or PETSC_NULL) 95470e92668SMatthew Knepley - ctx - location to stash Jacobian ctx (or PETSC_NULL) 955fee21e36SBarry Smith 95636851e7fSLois Curfman McInnes Level: advanced 95736851e7fSLois Curfman McInnes 958b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 959b4fd4287SBarry Smith @*/ 96070e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 961b4fd4287SBarry Smith { 9623a40ed3dSBarry Smith PetscFunctionBegin; 9634482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 964b4fd4287SBarry Smith if (A) *A = snes->jacobian; 965b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 96600036973SBarry Smith if (func) *func = snes->computejacobian; 96770e92668SMatthew Knepley if (ctx) *ctx = snes->jacP; 9683a40ed3dSBarry Smith PetscFunctionReturn(0); 969b4fd4287SBarry Smith } 970b4fd4287SBarry Smith 9719b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 97263dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 9739b94acceSBarry Smith 9744a2ae208SSatish Balay #undef __FUNCT__ 9754a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp" 9769b94acceSBarry Smith /*@ 9779b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 978272ac6f2SLois Curfman McInnes of a nonlinear solver. 9799b94acceSBarry Smith 980fee21e36SBarry Smith Collective on SNES 981fee21e36SBarry Smith 982c7afd0dbSLois Curfman McInnes Input Parameters: 98370e92668SMatthew Knepley . snes - the SNES context 984c7afd0dbSLois Curfman McInnes 985272ac6f2SLois Curfman McInnes Notes: 986272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 987272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 988272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 989272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 990272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 991272ac6f2SLois Curfman McInnes 99236851e7fSLois Curfman McInnes Level: advanced 99336851e7fSLois Curfman McInnes 9949b94acceSBarry Smith .keywords: SNES, nonlinear, setup 9959b94acceSBarry Smith 9969b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 9979b94acceSBarry Smith @*/ 99870e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 9999b94acceSBarry Smith { 1000dfbe8321SBarry Smith PetscErrorCode ierr; 10014b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 10023a40ed3dSBarry Smith 10033a40ed3dSBarry Smith PetscFunctionBegin; 10044482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 10054dc4c822SBarry Smith if (snes->setupcalled) PetscFunctionReturn(0); 10069b94acceSBarry Smith 1007b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1008c1f60f51SBarry Smith /* 1009c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1010dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1011c1f60f51SBarry Smith */ 1012c1f60f51SBarry Smith if (flg) { 1013c1f60f51SBarry Smith Mat J; 10145a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10155a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 101663ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator routines\n"));CHKERRQ(ierr); 10173a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 10183a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 1019c1f60f51SBarry Smith } 102045fc7adcSBarry Smith 1021c8d321e0SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) 102245fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 102345fc7adcSBarry Smith if (flg) { 102445fc7adcSBarry Smith Mat J; 102545fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 102645fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 102745fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 102845fc7adcSBarry Smith } 102932a4b47aSMatthew Knepley #endif 103045fc7adcSBarry Smith 1031b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1032c1f60f51SBarry Smith /* 1033dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1034c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1035c1f60f51SBarry Smith */ 103631615d04SLois Curfman McInnes if (flg) { 1037272ac6f2SLois Curfman McInnes Mat J; 1038b5d62d44SBarry Smith KSP ksp; 103994b7f48cSBarry Smith PC pc; 1040f3ef73ceSBarry Smith 10415a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10423a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 104363ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"));CHKERRQ(ierr); 10443a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 10453a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 10463a7fca6bSBarry Smith 1047f3ef73ceSBarry Smith /* force no preconditioner */ 104894b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1049b5d62d44SBarry Smith ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1050a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1051a9815358SBarry Smith if (!flg) { 1052f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1053272ac6f2SLois Curfman McInnes } 1054a9815358SBarry Smith } 1055f3ef73ceSBarry Smith 10561096aae1SMatthew Knepley if (!snes->vec_func && !snes->afine) { 10571096aae1SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 10581096aae1SMatthew Knepley } 10591096aae1SMatthew Knepley if (!snes->computefunction && !snes->afine) { 10601096aae1SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 10611096aae1SMatthew Knepley } 106229bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1063a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 106429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1065a8c6a408SBarry Smith } 1066a703fe33SLois Curfman McInnes 1067a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 10684b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 10696831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 107094b7f48cSBarry Smith KSP ksp; 107194b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1072186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1073a703fe33SLois Curfman McInnes } 10744b27c08aSLois Curfman McInnes 1075a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 10767aaed0d8SBarry Smith snes->setupcalled = PETSC_TRUE; 10773a40ed3dSBarry Smith PetscFunctionReturn(0); 10789b94acceSBarry Smith } 10799b94acceSBarry Smith 10804a2ae208SSatish Balay #undef __FUNCT__ 10814a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 10829b94acceSBarry Smith /*@C 10839b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 10849b94acceSBarry Smith with SNESCreate(). 10859b94acceSBarry Smith 1086c7afd0dbSLois Curfman McInnes Collective on SNES 1087c7afd0dbSLois Curfman McInnes 10889b94acceSBarry Smith Input Parameter: 10899b94acceSBarry Smith . snes - the SNES context 10909b94acceSBarry Smith 109136851e7fSLois Curfman McInnes Level: beginner 109236851e7fSLois Curfman McInnes 10939b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 10949b94acceSBarry Smith 109563a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 10969b94acceSBarry Smith @*/ 109763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 10989b94acceSBarry Smith { 109977431f27SBarry Smith PetscInt i; 11006849ba73SBarry Smith PetscErrorCode ierr; 11013a40ed3dSBarry Smith 11023a40ed3dSBarry Smith PetscFunctionBegin; 11034482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 11043a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1105d4bb536fSBarry Smith 1106be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 11070f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1108be0abb6dSBarry Smith 1109e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1110606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 11113a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 11123a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 11133ab0aad5SBarry Smith if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 111494b7f48cSBarry Smith ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1115522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1116b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1117b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1118b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1119b8a78c4aSBarry Smith } 1120b8a78c4aSBarry Smith } 1121a79aaaedSSatish Balay ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 11223a40ed3dSBarry Smith PetscFunctionReturn(0); 11239b94acceSBarry Smith } 11249b94acceSBarry Smith 11259b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 11269b94acceSBarry Smith 11274a2ae208SSatish Balay #undef __FUNCT__ 11284a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 11299b94acceSBarry Smith /*@ 1130d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 11319b94acceSBarry Smith 1132c7afd0dbSLois Curfman McInnes Collective on SNES 1133c7afd0dbSLois Curfman McInnes 11349b94acceSBarry Smith Input Parameters: 1135c7afd0dbSLois Curfman McInnes + snes - the SNES context 113670441072SBarry Smith . abstol - absolute convergence tolerance 113733174efeSLois Curfman McInnes . rtol - relative convergence tolerance 113833174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 113933174efeSLois Curfman McInnes of the change in the solution between steps 114033174efeSLois Curfman McInnes . maxit - maximum number of iterations 1141c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1142fee21e36SBarry Smith 114333174efeSLois Curfman McInnes Options Database Keys: 114470441072SBarry Smith + -snes_atol <abstol> - Sets abstol 1145c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1146c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1147c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1148c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 11499b94acceSBarry Smith 1150d7a720efSLois Curfman McInnes Notes: 11519b94acceSBarry Smith The default maximum number of iterations is 50. 11529b94acceSBarry Smith The default maximum number of function evaluations is 1000. 11539b94acceSBarry Smith 115436851e7fSLois Curfman McInnes Level: intermediate 115536851e7fSLois Curfman McInnes 115633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 11579b94acceSBarry Smith 11582492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance() 11599b94acceSBarry Smith @*/ 116063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 11619b94acceSBarry Smith { 11623a40ed3dSBarry Smith PetscFunctionBegin; 11634482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 116470441072SBarry Smith if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1165d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1166d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1167d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1168d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11693a40ed3dSBarry Smith PetscFunctionReturn(0); 11709b94acceSBarry Smith } 11719b94acceSBarry Smith 11724a2ae208SSatish Balay #undef __FUNCT__ 11734a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11749b94acceSBarry Smith /*@ 117533174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 117633174efeSLois Curfman McInnes 1177c7afd0dbSLois Curfman McInnes Not Collective 1178c7afd0dbSLois Curfman McInnes 117933174efeSLois Curfman McInnes Input Parameters: 1180c7afd0dbSLois Curfman McInnes + snes - the SNES context 118170441072SBarry Smith . abstol - absolute convergence tolerance 118233174efeSLois Curfman McInnes . rtol - relative convergence tolerance 118333174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 118433174efeSLois Curfman McInnes of the change in the solution between steps 118533174efeSLois Curfman McInnes . maxit - maximum number of iterations 1186c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1187fee21e36SBarry Smith 118833174efeSLois Curfman McInnes Notes: 118933174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 119033174efeSLois Curfman McInnes 119136851e7fSLois Curfman McInnes Level: intermediate 119236851e7fSLois Curfman McInnes 119333174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 119433174efeSLois Curfman McInnes 119533174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 119633174efeSLois Curfman McInnes @*/ 119763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 119833174efeSLois Curfman McInnes { 11993a40ed3dSBarry Smith PetscFunctionBegin; 12004482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 120170441072SBarry Smith if (abstol) *abstol = snes->abstol; 120233174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 120333174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 120433174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 120533174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 12063a40ed3dSBarry Smith PetscFunctionReturn(0); 120733174efeSLois Curfman McInnes } 120833174efeSLois Curfman McInnes 12094a2ae208SSatish Balay #undef __FUNCT__ 12104a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 121133174efeSLois Curfman McInnes /*@ 12129b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 12139b94acceSBarry Smith 1214fee21e36SBarry Smith Collective on SNES 1215fee21e36SBarry Smith 1216c7afd0dbSLois Curfman McInnes Input Parameters: 1217c7afd0dbSLois Curfman McInnes + snes - the SNES context 1218c7afd0dbSLois Curfman McInnes - tol - tolerance 1219c7afd0dbSLois Curfman McInnes 12209b94acceSBarry Smith Options Database Key: 1221c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 12229b94acceSBarry Smith 122336851e7fSLois Curfman McInnes Level: intermediate 122436851e7fSLois Curfman McInnes 12259b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 12269b94acceSBarry Smith 12272492ecdbSBarry Smith .seealso: SNESSetTolerances() 12289b94acceSBarry Smith @*/ 122963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 12309b94acceSBarry Smith { 12313a40ed3dSBarry Smith PetscFunctionBegin; 12324482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 12339b94acceSBarry Smith snes->deltatol = tol; 12343a40ed3dSBarry Smith PetscFunctionReturn(0); 12359b94acceSBarry Smith } 12369b94acceSBarry Smith 1237df9fa365SBarry Smith /* 1238df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1239df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1240df9fa365SBarry Smith macros instead of functions 1241df9fa365SBarry Smith */ 12424a2ae208SSatish Balay #undef __FUNCT__ 12434a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 124463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1245ce1608b8SBarry Smith { 1246dfbe8321SBarry Smith PetscErrorCode ierr; 1247ce1608b8SBarry Smith 1248ce1608b8SBarry Smith PetscFunctionBegin; 12494482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1250ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1251ce1608b8SBarry Smith PetscFunctionReturn(0); 1252ce1608b8SBarry Smith } 1253ce1608b8SBarry Smith 12544a2ae208SSatish Balay #undef __FUNCT__ 12554a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 125663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1257df9fa365SBarry Smith { 1258dfbe8321SBarry Smith PetscErrorCode ierr; 1259df9fa365SBarry Smith 1260df9fa365SBarry Smith PetscFunctionBegin; 1261df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1262df9fa365SBarry Smith PetscFunctionReturn(0); 1263df9fa365SBarry Smith } 1264df9fa365SBarry Smith 12654a2ae208SSatish Balay #undef __FUNCT__ 12664a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 126763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw) 1268df9fa365SBarry Smith { 1269dfbe8321SBarry Smith PetscErrorCode ierr; 1270df9fa365SBarry Smith 1271df9fa365SBarry Smith PetscFunctionBegin; 1272df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1273df9fa365SBarry Smith PetscFunctionReturn(0); 1274df9fa365SBarry Smith } 1275df9fa365SBarry Smith 12769b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12779b94acceSBarry Smith 12784a2ae208SSatish Balay #undef __FUNCT__ 12794a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12809b94acceSBarry Smith /*@C 12815cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12829b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12839b94acceSBarry Smith progress. 12849b94acceSBarry Smith 1285fee21e36SBarry Smith Collective on SNES 1286fee21e36SBarry Smith 1287c7afd0dbSLois Curfman McInnes Input Parameters: 1288c7afd0dbSLois Curfman McInnes + snes - the SNES context 1289c7afd0dbSLois Curfman McInnes . func - monitoring routine 1290b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1291b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desitre) 1292b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1293b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 12949b94acceSBarry Smith 1295c7afd0dbSLois Curfman McInnes Calling sequence of func: 1296a7cc72afSBarry Smith $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1297c7afd0dbSLois Curfman McInnes 1298c7afd0dbSLois Curfman McInnes + snes - the SNES context 1299c7afd0dbSLois Curfman McInnes . its - iteration number 1300c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 130140a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 13029b94acceSBarry Smith 13039665c990SLois Curfman McInnes Options Database Keys: 1304c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1305c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1306c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1307c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1308c7afd0dbSLois Curfman McInnes been hardwired into a code by 1309c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1310c7afd0dbSLois Curfman McInnes does not cancel those set via 1311c7afd0dbSLois Curfman McInnes the options database. 13129665c990SLois Curfman McInnes 1313639f9d9dSBarry Smith Notes: 13146bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 13156bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 13166bc08f3fSLois Curfman McInnes order in which they were set. 1317639f9d9dSBarry Smith 131836851e7fSLois Curfman McInnes Level: intermediate 131936851e7fSLois Curfman McInnes 13209b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 13219b94acceSBarry Smith 13225cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 13239b94acceSBarry Smith @*/ 132463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 13259b94acceSBarry Smith { 13263a40ed3dSBarry Smith PetscFunctionBegin; 13274482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1328639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 132929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1330639f9d9dSBarry Smith } 1331639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1332b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1333639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 13343a40ed3dSBarry Smith PetscFunctionReturn(0); 13359b94acceSBarry Smith } 13369b94acceSBarry Smith 13374a2ae208SSatish Balay #undef __FUNCT__ 13384a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 13395cd90555SBarry Smith /*@C 13405cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 13415cd90555SBarry Smith 1342c7afd0dbSLois Curfman McInnes Collective on SNES 1343c7afd0dbSLois Curfman McInnes 13445cd90555SBarry Smith Input Parameters: 13455cd90555SBarry Smith . snes - the SNES context 13465cd90555SBarry Smith 13471a480d89SAdministrator Options Database Key: 1348c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1349c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1350c7afd0dbSLois Curfman McInnes set via the options database 13515cd90555SBarry Smith 13525cd90555SBarry Smith Notes: 13535cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 13545cd90555SBarry Smith 135536851e7fSLois Curfman McInnes Level: intermediate 135636851e7fSLois Curfman McInnes 13575cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 13585cd90555SBarry Smith 13595cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 13605cd90555SBarry Smith @*/ 136163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes) 13625cd90555SBarry Smith { 13635cd90555SBarry Smith PetscFunctionBegin; 13644482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13655cd90555SBarry Smith snes->numbermonitors = 0; 13665cd90555SBarry Smith PetscFunctionReturn(0); 13675cd90555SBarry Smith } 13685cd90555SBarry Smith 13694a2ae208SSatish Balay #undef __FUNCT__ 13704a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13719b94acceSBarry Smith /*@C 13729b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13739b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13749b94acceSBarry Smith 1375fee21e36SBarry Smith Collective on SNES 1376fee21e36SBarry Smith 1377c7afd0dbSLois Curfman McInnes Input Parameters: 1378c7afd0dbSLois Curfman McInnes + snes - the SNES context 1379c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1380c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1381c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 13829b94acceSBarry Smith 1383c7afd0dbSLois Curfman McInnes Calling sequence of func: 13846849ba73SBarry Smith $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1385c7afd0dbSLois Curfman McInnes 1386c7afd0dbSLois Curfman McInnes + snes - the SNES context 1387c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1388184914b5SBarry Smith . reason - reason for convergence/divergence 1389c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 13904b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 13914b27c08aSLois Curfman McInnes - f - 2-norm of function 13929b94acceSBarry Smith 139336851e7fSLois Curfman McInnes Level: advanced 139436851e7fSLois Curfman McInnes 13959b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13969b94acceSBarry Smith 13974b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 13989b94acceSBarry Smith @*/ 139963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 14009b94acceSBarry Smith { 14013a40ed3dSBarry Smith PetscFunctionBegin; 14024482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14039b94acceSBarry Smith (snes)->converged = func; 14049b94acceSBarry Smith (snes)->cnvP = cctx; 14053a40ed3dSBarry Smith PetscFunctionReturn(0); 14069b94acceSBarry Smith } 14079b94acceSBarry Smith 14084a2ae208SSatish Balay #undef __FUNCT__ 14094a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 1410184914b5SBarry Smith /*@C 1411184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1412184914b5SBarry Smith 1413184914b5SBarry Smith Not Collective 1414184914b5SBarry Smith 1415184914b5SBarry Smith Input Parameter: 1416184914b5SBarry Smith . snes - the SNES context 1417184914b5SBarry Smith 1418184914b5SBarry Smith Output Parameter: 1419e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1420184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1421184914b5SBarry Smith 1422184914b5SBarry Smith Level: intermediate 1423184914b5SBarry Smith 1424184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1425184914b5SBarry Smith 1426184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1427184914b5SBarry Smith 14284b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1429184914b5SBarry Smith @*/ 143063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1431184914b5SBarry Smith { 1432184914b5SBarry Smith PetscFunctionBegin; 14334482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14344482741eSBarry Smith PetscValidPointer(reason,2); 1435184914b5SBarry Smith *reason = snes->reason; 1436184914b5SBarry Smith PetscFunctionReturn(0); 1437184914b5SBarry Smith } 1438184914b5SBarry Smith 14394a2ae208SSatish Balay #undef __FUNCT__ 14404a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1441c9005455SLois Curfman McInnes /*@ 1442c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1443c9005455SLois Curfman McInnes 1444fee21e36SBarry Smith Collective on SNES 1445fee21e36SBarry Smith 1446c7afd0dbSLois Curfman McInnes Input Parameters: 1447c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1448c7afd0dbSLois Curfman McInnes . a - array to hold history 1449cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1450758f92a0SBarry Smith . na - size of a and its 145164731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1452758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1453c7afd0dbSLois Curfman McInnes 1454c9005455SLois Curfman McInnes Notes: 14554b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1456c9005455SLois Curfman McInnes at each step. 1457c9005455SLois Curfman McInnes 1458c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1459c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1460c9005455SLois Curfman McInnes during the section of code that is being timed. 1461c9005455SLois Curfman McInnes 146236851e7fSLois Curfman McInnes Level: intermediate 146336851e7fSLois Curfman McInnes 1464c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1465758f92a0SBarry Smith 146608405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1467758f92a0SBarry Smith 1468c9005455SLois Curfman McInnes @*/ 146963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1470c9005455SLois Curfman McInnes { 14713a40ed3dSBarry Smith PetscFunctionBegin; 14724482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14734482741eSBarry Smith if (na) PetscValidScalarPointer(a,2); 1474c9005455SLois Curfman McInnes snes->conv_hist = a; 1475758f92a0SBarry Smith snes->conv_hist_its = its; 1476758f92a0SBarry Smith snes->conv_hist_max = na; 1477758f92a0SBarry Smith snes->conv_hist_reset = reset; 1478758f92a0SBarry Smith PetscFunctionReturn(0); 1479758f92a0SBarry Smith } 1480758f92a0SBarry Smith 14814a2ae208SSatish Balay #undef __FUNCT__ 14824a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 14830c4c9dddSBarry Smith /*@C 1484758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1485758f92a0SBarry Smith 1486758f92a0SBarry Smith Collective on SNES 1487758f92a0SBarry Smith 1488758f92a0SBarry Smith Input Parameter: 1489758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1490758f92a0SBarry Smith 1491758f92a0SBarry Smith Output Parameters: 1492758f92a0SBarry Smith . a - array to hold history 1493758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1494758f92a0SBarry Smith negative if not converged) for each solve. 1495758f92a0SBarry Smith - na - size of a and its 1496758f92a0SBarry Smith 1497758f92a0SBarry Smith Notes: 1498758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1499758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1500758f92a0SBarry Smith 1501758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1502758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1503758f92a0SBarry Smith during the section of code that is being timed. 1504758f92a0SBarry Smith 1505758f92a0SBarry Smith Level: intermediate 1506758f92a0SBarry Smith 1507758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1508758f92a0SBarry Smith 1509758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1510758f92a0SBarry Smith 1511758f92a0SBarry Smith @*/ 151263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1513758f92a0SBarry Smith { 1514758f92a0SBarry Smith PetscFunctionBegin; 15154482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1516758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1517758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1518758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 15193a40ed3dSBarry Smith PetscFunctionReturn(0); 1520c9005455SLois Curfman McInnes } 1521c9005455SLois Curfman McInnes 1522e74ef692SMatthew Knepley #undef __FUNCT__ 1523e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 1524ac226902SBarry Smith /*@C 152576b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 152676b2cf59SMatthew Knepley at the beginning of every step of the iteration. 152776b2cf59SMatthew Knepley 152876b2cf59SMatthew Knepley Collective on SNES 152976b2cf59SMatthew Knepley 153076b2cf59SMatthew Knepley Input Parameters: 153176b2cf59SMatthew Knepley . snes - The nonlinear solver context 153276b2cf59SMatthew Knepley . func - The function 153376b2cf59SMatthew Knepley 153476b2cf59SMatthew Knepley Calling sequence of func: 1535b5d30489SBarry Smith . func (SNES snes, PetscInt step); 153676b2cf59SMatthew Knepley 153776b2cf59SMatthew Knepley . step - The current step of the iteration 153876b2cf59SMatthew Knepley 153976b2cf59SMatthew Knepley Level: intermediate 154076b2cf59SMatthew Knepley 154176b2cf59SMatthew Knepley .keywords: SNES, update 1542b5d30489SBarry Smith 154376b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 154476b2cf59SMatthew Knepley @*/ 154563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 154676b2cf59SMatthew Knepley { 154776b2cf59SMatthew Knepley PetscFunctionBegin; 15484482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 154976b2cf59SMatthew Knepley snes->update = func; 155076b2cf59SMatthew Knepley PetscFunctionReturn(0); 155176b2cf59SMatthew Knepley } 155276b2cf59SMatthew Knepley 1553e74ef692SMatthew Knepley #undef __FUNCT__ 1554e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 155576b2cf59SMatthew Knepley /*@ 155676b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 155776b2cf59SMatthew Knepley 155876b2cf59SMatthew Knepley Not collective 155976b2cf59SMatthew Knepley 156076b2cf59SMatthew Knepley Input Parameters: 156176b2cf59SMatthew Knepley . snes - The nonlinear solver context 156276b2cf59SMatthew Knepley . step - The current step of the iteration 156376b2cf59SMatthew Knepley 1564205452f4SMatthew Knepley Level: intermediate 1565205452f4SMatthew Knepley 156676b2cf59SMatthew Knepley .keywords: SNES, update 156776b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 156876b2cf59SMatthew Knepley @*/ 156963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 157076b2cf59SMatthew Knepley { 157176b2cf59SMatthew Knepley PetscFunctionBegin; 157276b2cf59SMatthew Knepley PetscFunctionReturn(0); 157376b2cf59SMatthew Knepley } 157476b2cf59SMatthew Knepley 15754a2ae208SSatish Balay #undef __FUNCT__ 15764a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 15779b94acceSBarry Smith /* 15789b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 15799b94acceSBarry Smith positive parameter delta. 15809b94acceSBarry Smith 15819b94acceSBarry Smith Input Parameters: 1582c7afd0dbSLois Curfman McInnes + snes - the SNES context 15839b94acceSBarry Smith . y - approximate solution of linear system 15849b94acceSBarry Smith . fnorm - 2-norm of current function 1585c7afd0dbSLois Curfman McInnes - delta - trust region size 15869b94acceSBarry Smith 15879b94acceSBarry Smith Output Parameters: 1588c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 15899b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 15909b94acceSBarry Smith region, and exceeds zero otherwise. 1591c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 15929b94acceSBarry Smith 15939b94acceSBarry Smith Note: 15944b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 15959b94acceSBarry Smith is set to be the maximum allowable step size. 15969b94acceSBarry Smith 15979b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 15989b94acceSBarry Smith */ 1599dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 16009b94acceSBarry Smith { 1601064f8208SBarry Smith PetscReal nrm; 1602ea709b57SSatish Balay PetscScalar cnorm; 1603dfbe8321SBarry Smith PetscErrorCode ierr; 16043a40ed3dSBarry Smith 16053a40ed3dSBarry Smith PetscFunctionBegin; 16064482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16074482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1608c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,2); 1609184914b5SBarry Smith 1610064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1611064f8208SBarry Smith if (nrm > *delta) { 1612064f8208SBarry Smith nrm = *delta/nrm; 1613064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1614064f8208SBarry Smith cnorm = nrm; 16152dcb1b2aSMatthew Knepley ierr = VecScale(y,cnorm);CHKERRQ(ierr); 16169b94acceSBarry Smith *ynorm = *delta; 16179b94acceSBarry Smith } else { 16189b94acceSBarry Smith *gpnorm = 0.0; 1619064f8208SBarry Smith *ynorm = nrm; 16209b94acceSBarry Smith } 16213a40ed3dSBarry Smith PetscFunctionReturn(0); 16229b94acceSBarry Smith } 16239b94acceSBarry Smith 16244a2ae208SSatish Balay #undef __FUNCT__ 16254a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 16269b94acceSBarry Smith /*@ 1627f69a0ea3SMatthew Knepley SNESSolve - Solves a nonlinear system F(x) = b. 1628f69a0ea3SMatthew Knepley Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 16299b94acceSBarry Smith 1630c7afd0dbSLois Curfman McInnes Collective on SNES 1631c7afd0dbSLois Curfman McInnes 1632b2002411SLois Curfman McInnes Input Parameters: 1633c7afd0dbSLois Curfman McInnes + snes - the SNES context 1634f69a0ea3SMatthew Knepley . b - the constant part of the equation, or PETSC_NULL to use zero. 163570e92668SMatthew Knepley - x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution() 16369b94acceSBarry Smith 1637b2002411SLois Curfman McInnes Notes: 16388ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 16398ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 16408ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 16418ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 16428ddd3da0SLois Curfman McInnes 164336851e7fSLois Curfman McInnes Level: beginner 164436851e7fSLois Curfman McInnes 16459b94acceSBarry Smith .keywords: SNES, nonlinear, solve 16469b94acceSBarry Smith 164770e92668SMatthew Knepley .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution() 16489b94acceSBarry Smith @*/ 1649f69a0ea3SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 16509b94acceSBarry Smith { 1651dfbe8321SBarry Smith PetscErrorCode ierr; 1652f1af5d2fSBarry Smith PetscTruth flg; 1653052efed2SBarry Smith 16543a40ed3dSBarry Smith PetscFunctionBegin; 16554482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16561302d50aSBarry Smith if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 165774637425SBarry Smith 1658f69a0ea3SMatthew Knepley if (b) { 1659f69a0ea3SMatthew Knepley ierr = SNESSetRhs(snes, b); CHKERRQ(ierr); 16601096aae1SMatthew Knepley if (!snes->vec_func) { 16611096aae1SMatthew Knepley Vec r; 16621096aae1SMatthew Knepley 16631096aae1SMatthew Knepley ierr = VecDuplicate(b, &r); CHKERRQ(ierr); 16641096aae1SMatthew Knepley ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr); 16651096aae1SMatthew Knepley } 1666f69a0ea3SMatthew Knepley } 166770e92668SMatthew Knepley if (x) { 1668f69a0ea3SMatthew Knepley PetscValidHeaderSpecific(x,VEC_COOKIE,3); 1669f69a0ea3SMatthew Knepley PetscCheckSameComm(snes,1,x,3); 167070e92668SMatthew Knepley } else { 167170e92668SMatthew Knepley ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr); 167270e92668SMatthew Knepley if (!x) { 167370e92668SMatthew Knepley ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr); 167470e92668SMatthew Knepley } 167570e92668SMatthew Knepley } 167670e92668SMatthew Knepley snes->vec_sol = snes->vec_sol_always = x; 167770e92668SMatthew Knepley if (!snes->setupcalled) { 167870e92668SMatthew Knepley ierr = SNESSetUp(snes);CHKERRQ(ierr); 167970e92668SMatthew Knepley } 1680abc0a331SBarry Smith if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1681d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 168250ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1683d5e45103SBarry Smith 1684d5e45103SBarry Smith ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN); 1685d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 168638f152feSBarry Smith /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */ 168719717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr); 1688d5e45103SBarry Smith } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) { 1689d5e45103SBarry Smith /* translate exception into SNES not converged reason */ 1690d5e45103SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 169138f152feSBarry Smith ierr = 0; 169238f152feSBarry Smith } 169338f152feSBarry Smith CHKERRQ(ierr); 1694d5e45103SBarry Smith 1695d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1696b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 16978b179ff8SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1698da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1699da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 17005968eb51SBarry Smith if (snes->printreason) { 17015968eb51SBarry Smith if (snes->reason > 0) { 17029dcbbd2bSBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 17035968eb51SBarry Smith } else { 17049dcbbd2bSBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 17055968eb51SBarry Smith } 17065968eb51SBarry Smith } 17075968eb51SBarry Smith 17083a40ed3dSBarry Smith PetscFunctionReturn(0); 17099b94acceSBarry Smith } 17109b94acceSBarry Smith 17119b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 17129b94acceSBarry Smith 17134a2ae208SSatish Balay #undef __FUNCT__ 17144a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 171582bf6240SBarry Smith /*@C 17164b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 17179b94acceSBarry Smith 1718fee21e36SBarry Smith Collective on SNES 1719fee21e36SBarry Smith 1720c7afd0dbSLois Curfman McInnes Input Parameters: 1721c7afd0dbSLois Curfman McInnes + snes - the SNES context 1722454a90a3SBarry Smith - type - a known method 1723c7afd0dbSLois Curfman McInnes 1724c7afd0dbSLois Curfman McInnes Options Database Key: 1725454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1726c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1727ae12b187SLois Curfman McInnes 17289b94acceSBarry Smith Notes: 1729e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 17304b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1731c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17324b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1733c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17349b94acceSBarry Smith 1735ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1736ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1737ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1738ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1739ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1740ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1741ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1742ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1743ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1744b0a32e0cSBarry Smith appropriate method. 174536851e7fSLois Curfman McInnes 174636851e7fSLois Curfman McInnes Level: intermediate 1747a703fe33SLois Curfman McInnes 1748454a90a3SBarry Smith .keywords: SNES, set, type 1749435da068SBarry Smith 1750435da068SBarry Smith .seealso: SNESType, SNESCreate() 1751435da068SBarry Smith 17529b94acceSBarry Smith @*/ 1753e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type) 17549b94acceSBarry Smith { 1755dfbe8321SBarry Smith PetscErrorCode ierr,(*r)(SNES); 17566831982aSBarry Smith PetscTruth match; 17573a40ed3dSBarry Smith 17583a40ed3dSBarry Smith PetscFunctionBegin; 17594482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17604482741eSBarry Smith PetscValidCharPointer(type,2); 176182bf6240SBarry Smith 17626831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 17630f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 176482bf6240SBarry Smith 176582bf6240SBarry Smith if (snes->setupcalled) { 17664dc4c822SBarry Smith snes->setupcalled = PETSC_FALSE; 1767e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 176882bf6240SBarry Smith snes->data = 0; 176982bf6240SBarry Smith } 17707d1a2b2bSBarry Smith 17719b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 177282bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1773b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1774958c9bccSBarry Smith if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1775606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 177682bf6240SBarry Smith snes->data = 0; 17773a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1778454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 17793a40ed3dSBarry Smith PetscFunctionReturn(0); 17809b94acceSBarry Smith } 17819b94acceSBarry Smith 1782a847f771SSatish Balay 17839b94acceSBarry Smith /* --------------------------------------------------------------------- */ 17844a2ae208SSatish Balay #undef __FUNCT__ 17854a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 17869b94acceSBarry Smith /*@C 17879b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1788f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 17899b94acceSBarry Smith 1790fee21e36SBarry Smith Not Collective 1791fee21e36SBarry Smith 179236851e7fSLois Curfman McInnes Level: advanced 179336851e7fSLois Curfman McInnes 17949b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 17959b94acceSBarry Smith 17969b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 17979b94acceSBarry Smith @*/ 179863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 17999b94acceSBarry Smith { 1800dfbe8321SBarry Smith PetscErrorCode ierr; 180182bf6240SBarry Smith 18023a40ed3dSBarry Smith PetscFunctionBegin; 180382bf6240SBarry Smith if (SNESList) { 1804b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 180582bf6240SBarry Smith SNESList = 0; 18069b94acceSBarry Smith } 18074c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 18083a40ed3dSBarry Smith PetscFunctionReturn(0); 18099b94acceSBarry Smith } 18109b94acceSBarry Smith 18114a2ae208SSatish Balay #undef __FUNCT__ 18124a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 18139b94acceSBarry Smith /*@C 18149a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 18159b94acceSBarry Smith 1816c7afd0dbSLois Curfman McInnes Not Collective 1817c7afd0dbSLois Curfman McInnes 18189b94acceSBarry Smith Input Parameter: 18194b0e389bSBarry Smith . snes - nonlinear solver context 18209b94acceSBarry Smith 18219b94acceSBarry Smith Output Parameter: 18223a7fca6bSBarry Smith . type - SNES method (a character string) 18239b94acceSBarry Smith 182436851e7fSLois Curfman McInnes Level: intermediate 182536851e7fSLois Curfman McInnes 1826454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 18279b94acceSBarry Smith @*/ 182863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type) 18299b94acceSBarry Smith { 18303a40ed3dSBarry Smith PetscFunctionBegin; 18314482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18324482741eSBarry Smith PetscValidPointer(type,2); 1833454a90a3SBarry Smith *type = snes->type_name; 18343a40ed3dSBarry Smith PetscFunctionReturn(0); 18359b94acceSBarry Smith } 18369b94acceSBarry Smith 18374a2ae208SSatish Balay #undef __FUNCT__ 18384a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 18399b94acceSBarry Smith /*@C 18409b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 18419b94acceSBarry Smith stored. 18429b94acceSBarry Smith 1843c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1844c7afd0dbSLois Curfman McInnes 18459b94acceSBarry Smith Input Parameter: 18469b94acceSBarry Smith . snes - the SNES context 18479b94acceSBarry Smith 18489b94acceSBarry Smith Output Parameter: 18499b94acceSBarry Smith . x - the solution 18509b94acceSBarry Smith 185170e92668SMatthew Knepley Level: intermediate 185236851e7fSLois Curfman McInnes 18539b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 18549b94acceSBarry Smith 185570e92668SMatthew Knepley .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 18569b94acceSBarry Smith @*/ 185763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 18589b94acceSBarry Smith { 18593a40ed3dSBarry Smith PetscFunctionBegin; 18604482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18614482741eSBarry Smith PetscValidPointer(x,2); 18629b94acceSBarry Smith *x = snes->vec_sol_always; 18633a40ed3dSBarry Smith PetscFunctionReturn(0); 18649b94acceSBarry Smith } 18659b94acceSBarry Smith 18664a2ae208SSatish Balay #undef __FUNCT__ 186770e92668SMatthew Knepley #define __FUNCT__ "SNESSetSolution" 186870e92668SMatthew Knepley /*@C 186970e92668SMatthew Knepley SNESSetSolution - Sets the vector where the approximate solution is stored. 187070e92668SMatthew Knepley 187170e92668SMatthew Knepley Not Collective, but Vec is parallel if SNES is parallel 187270e92668SMatthew Knepley 187370e92668SMatthew Knepley Input Parameters: 187470e92668SMatthew Knepley + snes - the SNES context 187570e92668SMatthew Knepley - x - the solution 187670e92668SMatthew Knepley 187770e92668SMatthew Knepley Output Parameter: 187870e92668SMatthew Knepley 187970e92668SMatthew Knepley Level: intermediate 188070e92668SMatthew Knepley 188170e92668SMatthew Knepley .keywords: SNES, nonlinear, set, solution 188270e92668SMatthew Knepley 188370e92668SMatthew Knepley .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 188470e92668SMatthew Knepley @*/ 188570e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x) 188670e92668SMatthew Knepley { 188770e92668SMatthew Knepley PetscFunctionBegin; 188870e92668SMatthew Knepley PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 188970e92668SMatthew Knepley PetscValidHeaderSpecific(x,VEC_COOKIE,2); 189070e92668SMatthew Knepley PetscCheckSameComm(snes,1,x,2); 189170e92668SMatthew Knepley snes->vec_sol_always = x; 189270e92668SMatthew Knepley PetscFunctionReturn(0); 189370e92668SMatthew Knepley } 189470e92668SMatthew Knepley 189570e92668SMatthew Knepley #undef __FUNCT__ 18964a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 18979b94acceSBarry Smith /*@C 18989b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 18999b94acceSBarry Smith stored. 19009b94acceSBarry Smith 1901c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1902c7afd0dbSLois Curfman McInnes 19039b94acceSBarry Smith Input Parameter: 19049b94acceSBarry Smith . snes - the SNES context 19059b94acceSBarry Smith 19069b94acceSBarry Smith Output Parameter: 19079b94acceSBarry Smith . x - the solution update 19089b94acceSBarry Smith 190936851e7fSLois Curfman McInnes Level: advanced 191036851e7fSLois Curfman McInnes 19119b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 19129b94acceSBarry Smith 19139b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 19149b94acceSBarry Smith @*/ 191563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 19169b94acceSBarry Smith { 19173a40ed3dSBarry Smith PetscFunctionBegin; 19184482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19194482741eSBarry Smith PetscValidPointer(x,2); 19209b94acceSBarry Smith *x = snes->vec_sol_update_always; 19213a40ed3dSBarry Smith PetscFunctionReturn(0); 19229b94acceSBarry Smith } 19239b94acceSBarry Smith 19244a2ae208SSatish Balay #undef __FUNCT__ 19254a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 19269b94acceSBarry Smith /*@C 19273638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 19289b94acceSBarry Smith 1929c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1930c7afd0dbSLois Curfman McInnes 19319b94acceSBarry Smith Input Parameter: 19329b94acceSBarry Smith . snes - the SNES context 19339b94acceSBarry Smith 19349b94acceSBarry Smith Output Parameter: 19357bf4e008SBarry Smith + r - the function (or PETSC_NULL) 193670e92668SMatthew Knepley . func - the function (or PETSC_NULL) 193770e92668SMatthew Knepley - ctx - the function context (or PETSC_NULL) 19389b94acceSBarry Smith 193936851e7fSLois Curfman McInnes Level: advanced 194036851e7fSLois Curfman McInnes 1941a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 19429b94acceSBarry Smith 19434b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 19449b94acceSBarry Smith @*/ 194570e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 19469b94acceSBarry Smith { 19473a40ed3dSBarry Smith PetscFunctionBegin; 19484482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19497bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 195000036973SBarry Smith if (func) *func = snes->computefunction; 195170e92668SMatthew Knepley if (ctx) *ctx = snes->funP; 19523a40ed3dSBarry Smith PetscFunctionReturn(0); 19539b94acceSBarry Smith } 19549b94acceSBarry Smith 19554a2ae208SSatish Balay #undef __FUNCT__ 19564a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 19573c7409f5SSatish Balay /*@C 19583c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1959d850072dSLois Curfman McInnes SNES options in the database. 19603c7409f5SSatish Balay 1961fee21e36SBarry Smith Collective on SNES 1962fee21e36SBarry Smith 1963c7afd0dbSLois Curfman McInnes Input Parameter: 1964c7afd0dbSLois Curfman McInnes + snes - the SNES context 1965c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1966c7afd0dbSLois Curfman McInnes 1967d850072dSLois Curfman McInnes Notes: 1968a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1969c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1970d850072dSLois Curfman McInnes 197136851e7fSLois Curfman McInnes Level: advanced 197236851e7fSLois Curfman McInnes 19733c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1974a86d99e1SLois Curfman McInnes 1975a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19763c7409f5SSatish Balay @*/ 197763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 19783c7409f5SSatish Balay { 1979dfbe8321SBarry Smith PetscErrorCode ierr; 19803c7409f5SSatish Balay 19813a40ed3dSBarry Smith PetscFunctionBegin; 19824482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1983639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 198494b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 19853a40ed3dSBarry Smith PetscFunctionReturn(0); 19863c7409f5SSatish Balay } 19873c7409f5SSatish Balay 19884a2ae208SSatish Balay #undef __FUNCT__ 19894a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 19903c7409f5SSatish Balay /*@C 1991f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1992d850072dSLois Curfman McInnes SNES options in the database. 19933c7409f5SSatish Balay 1994fee21e36SBarry Smith Collective on SNES 1995fee21e36SBarry Smith 1996c7afd0dbSLois Curfman McInnes Input Parameters: 1997c7afd0dbSLois Curfman McInnes + snes - the SNES context 1998c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1999c7afd0dbSLois Curfman McInnes 2000d850072dSLois Curfman McInnes Notes: 2001a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2002c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2003d850072dSLois Curfman McInnes 200436851e7fSLois Curfman McInnes Level: advanced 200536851e7fSLois Curfman McInnes 20063c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2007a86d99e1SLois Curfman McInnes 2008a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20093c7409f5SSatish Balay @*/ 201063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 20113c7409f5SSatish Balay { 2012dfbe8321SBarry Smith PetscErrorCode ierr; 20133c7409f5SSatish Balay 20143a40ed3dSBarry Smith PetscFunctionBegin; 20154482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2016639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 201794b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20183a40ed3dSBarry Smith PetscFunctionReturn(0); 20193c7409f5SSatish Balay } 20203c7409f5SSatish Balay 20214a2ae208SSatish Balay #undef __FUNCT__ 20224a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 20239ab63eb5SSatish Balay /*@C 20243c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20253c7409f5SSatish Balay SNES options in the database. 20263c7409f5SSatish Balay 2027c7afd0dbSLois Curfman McInnes Not Collective 2028c7afd0dbSLois Curfman McInnes 20293c7409f5SSatish Balay Input Parameter: 20303c7409f5SSatish Balay . snes - the SNES context 20313c7409f5SSatish Balay 20323c7409f5SSatish Balay Output Parameter: 20333c7409f5SSatish Balay . prefix - pointer to the prefix string used 20343c7409f5SSatish Balay 20359ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 20369ab63eb5SSatish Balay sufficient length to hold the prefix. 20379ab63eb5SSatish Balay 203836851e7fSLois Curfman McInnes Level: advanced 203936851e7fSLois Curfman McInnes 20403c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2041a86d99e1SLois Curfman McInnes 2042a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20433c7409f5SSatish Balay @*/ 2044e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 20453c7409f5SSatish Balay { 2046dfbe8321SBarry Smith PetscErrorCode ierr; 20473c7409f5SSatish Balay 20483a40ed3dSBarry Smith PetscFunctionBegin; 20494482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2050639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 20513a40ed3dSBarry Smith PetscFunctionReturn(0); 20523c7409f5SSatish Balay } 20533c7409f5SSatish Balay 2054b2002411SLois Curfman McInnes 20554a2ae208SSatish Balay #undef __FUNCT__ 20564a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 20573cea93caSBarry Smith /*@C 20583cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 20593cea93caSBarry Smith 20607f6c08e0SMatthew Knepley Level: advanced 20613cea93caSBarry Smith @*/ 206263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2063b2002411SLois Curfman McInnes { 2064e2d1d2b7SBarry Smith char fullname[PETSC_MAX_PATH_LEN]; 2065dfbe8321SBarry Smith PetscErrorCode ierr; 2066b2002411SLois Curfman McInnes 2067b2002411SLois Curfman McInnes PetscFunctionBegin; 2068b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2069c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2070b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2071b2002411SLois Curfman McInnes } 2072da9b6338SBarry Smith 2073da9b6338SBarry Smith #undef __FUNCT__ 2074da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 207563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2076da9b6338SBarry Smith { 2077dfbe8321SBarry Smith PetscErrorCode ierr; 207877431f27SBarry Smith PetscInt N,i,j; 2079da9b6338SBarry Smith Vec u,uh,fh; 2080da9b6338SBarry Smith PetscScalar value; 2081da9b6338SBarry Smith PetscReal norm; 2082da9b6338SBarry Smith 2083da9b6338SBarry Smith PetscFunctionBegin; 2084da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2085da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2086da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2087da9b6338SBarry Smith 2088da9b6338SBarry Smith /* currently only works for sequential */ 2089da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2090da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2091da9b6338SBarry Smith for (i=0; i<N; i++) { 2092da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 209377431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2094da9b6338SBarry Smith for (j=-10; j<11; j++) { 2095ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2096da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 20973ab0aad5SBarry Smith ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2098da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 209977431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2100da9b6338SBarry Smith value = -value; 2101da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2102da9b6338SBarry Smith } 2103da9b6338SBarry Smith } 2104da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2105da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2106da9b6338SBarry Smith PetscFunctionReturn(0); 2107da9b6338SBarry Smith } 2108