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; 48454a90a3SBarry Smith char *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 225b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr); 226186905e3SBarry Smith 227b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 22887828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 22987828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 23087828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 23187828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 23287828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 23387828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 234186905e3SBarry Smith 235b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr); 23693c39befSBarry Smith if (flg) {snes->converged = 0;} 237b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr); 2385cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 239b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr); 240b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);} 2413a7fca6bSBarry Smith ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr); 2423a7fca6bSBarry Smith if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);} 243b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr); 244b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);} 245b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr); 246b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 247b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr); 2487c922b88SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);} 2495ed2d596SBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr); 2505ed2d596SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);} 251b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr); 252186905e3SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 253e24b481bSBarry Smith 254b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 2554b27c08aSLois Curfman McInnes if (flg) { 256186905e3SBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 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; 5884482741eSBarry Smith PetscValidPointer(outsnes,1); 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; 61182bf6240SBarry Smith snes->setupcalled = 0; 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 7093ab0aad5SBarry Smith .seealso: 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__ 7304a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction" 7319b94acceSBarry Smith /*@ 73236851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7339b94acceSBarry Smith SNESSetFunction(). 7349b94acceSBarry Smith 735c7afd0dbSLois Curfman McInnes Collective on SNES 736c7afd0dbSLois Curfman McInnes 7379b94acceSBarry Smith Input Parameters: 738c7afd0dbSLois Curfman McInnes + snes - the SNES context 739c7afd0dbSLois Curfman McInnes - x - input vector 7409b94acceSBarry Smith 7419b94acceSBarry Smith Output Parameter: 7423638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7439b94acceSBarry Smith 7441bffabb2SLois Curfman McInnes Notes: 74536851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 74636851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 74736851e7fSLois Curfman McInnes themselves. 74836851e7fSLois Curfman McInnes 74936851e7fSLois Curfman McInnes Level: developer 75036851e7fSLois Curfman McInnes 7519b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 7529b94acceSBarry Smith 753a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 7549b94acceSBarry Smith @*/ 75563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y) 7569b94acceSBarry Smith { 757dfbe8321SBarry Smith PetscErrorCode ierr; 7589b94acceSBarry Smith 7593a40ed3dSBarry Smith PetscFunctionBegin; 7604482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7614482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 7624482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,3); 763c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 764c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,3); 765184914b5SBarry Smith 766d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 767d64ed03dSBarry Smith PetscStackPush("SNES user function"); 76819717074SBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); 769d64ed03dSBarry Smith PetscStackPop; 770*d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 77119717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 77219717074SBarry Smith } 773*d5e45103SBarry Smith CHKERRQ(ierr); 7743ab0aad5SBarry Smith if (snes->afine) { 7753ab0aad5SBarry Smith PetscScalar mone = -1.0; 7763ab0aad5SBarry Smith ierr = VecAXPY(&mone,snes->afine,y);CHKERRQ(ierr); 7773ab0aad5SBarry Smith } 778ae3c334cSLois Curfman McInnes snes->nfuncs++; 779d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 7803a40ed3dSBarry Smith PetscFunctionReturn(0); 7819b94acceSBarry Smith } 7829b94acceSBarry Smith 7834a2ae208SSatish Balay #undef __FUNCT__ 7844a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian" 78562fef451SLois Curfman McInnes /*@ 78662fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 78762fef451SLois Curfman McInnes set with SNESSetJacobian(). 78862fef451SLois Curfman McInnes 789c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 790c7afd0dbSLois Curfman McInnes 79162fef451SLois Curfman McInnes Input Parameters: 792c7afd0dbSLois Curfman McInnes + snes - the SNES context 793c7afd0dbSLois Curfman McInnes - x - input vector 79462fef451SLois Curfman McInnes 79562fef451SLois Curfman McInnes Output Parameters: 796c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 79762fef451SLois Curfman McInnes . B - optional preconditioning matrix 798c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 799fee21e36SBarry Smith 80062fef451SLois Curfman McInnes Notes: 80162fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 80262fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 80362fef451SLois Curfman McInnes 80494b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 805dc5a77f8SLois Curfman McInnes flag parameter. 80662fef451SLois Curfman McInnes 80736851e7fSLois Curfman McInnes Level: developer 80836851e7fSLois Curfman McInnes 80962fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 81062fef451SLois Curfman McInnes 81194b7f48cSBarry Smith .seealso: SNESSetJacobian(), KSPSetOperators() 81262fef451SLois Curfman McInnes @*/ 81363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8149b94acceSBarry Smith { 815dfbe8321SBarry Smith PetscErrorCode ierr; 8163a40ed3dSBarry Smith 8173a40ed3dSBarry Smith PetscFunctionBegin; 8184482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8194482741eSBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE,2); 8204482741eSBarry Smith PetscValidPointer(flg,5); 821c9780b6fSBarry Smith PetscCheckSameComm(snes,1,X,2); 8223a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 823d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 824c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 825d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 8269b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 827d64ed03dSBarry Smith PetscStackPop; 828d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 8296d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 8304482741eSBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 8314482741eSBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 8323a40ed3dSBarry Smith PetscFunctionReturn(0); 8339b94acceSBarry Smith } 8349b94acceSBarry Smith 8354a2ae208SSatish Balay #undef __FUNCT__ 8364a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian" 8379b94acceSBarry Smith /*@C 8389b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 839044dda88SLois Curfman McInnes location to store the matrix. 8409b94acceSBarry Smith 841c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 842c7afd0dbSLois Curfman McInnes 8439b94acceSBarry Smith Input Parameters: 844c7afd0dbSLois Curfman McInnes + snes - the SNES context 8459b94acceSBarry Smith . A - Jacobian matrix 8469b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 8479b94acceSBarry Smith . func - Jacobian evaluation routine 848c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 8492cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 8509b94acceSBarry Smith 8519b94acceSBarry Smith Calling sequence of func: 8528d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 8539b94acceSBarry Smith 854c7afd0dbSLois Curfman McInnes + x - input vector 8559b94acceSBarry Smith . A - Jacobian matrix 8569b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 857ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 85894b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 859c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 8609b94acceSBarry Smith 8619b94acceSBarry Smith Notes: 86294b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 8632cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 864ac21db08SLois Curfman McInnes 865ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 8669b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 8679b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 8689b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 8699b94acceSBarry Smith throughout the global iterations. 8709b94acceSBarry Smith 87136851e7fSLois Curfman McInnes Level: beginner 87236851e7fSLois Curfman McInnes 8739b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 8749b94acceSBarry Smith 87513f74950SBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), , MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor() 8769b94acceSBarry Smith @*/ 87763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 8789b94acceSBarry Smith { 879dfbe8321SBarry Smith PetscErrorCode ierr; 8803a7fca6bSBarry Smith 8813a40ed3dSBarry Smith PetscFunctionBegin; 8824482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8834482741eSBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 8844482741eSBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 885c9780b6fSBarry Smith if (A) PetscCheckSameComm(snes,1,A,2); 886c9780b6fSBarry Smith if (B) PetscCheckSameComm(snes,1,B,2); 8873a7fca6bSBarry Smith if (func) snes->computejacobian = func; 8883a7fca6bSBarry Smith if (ctx) snes->jacP = ctx; 8893a7fca6bSBarry Smith if (A) { 8903a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 8919b94acceSBarry Smith snes->jacobian = A; 8923a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 8933a7fca6bSBarry Smith } 8943a7fca6bSBarry Smith if (B) { 8953a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 8969b94acceSBarry Smith snes->jacobian_pre = B; 8973a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 8983a7fca6bSBarry Smith } 89969a612a9SBarry Smith ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 9003a40ed3dSBarry Smith PetscFunctionReturn(0); 9019b94acceSBarry Smith } 90262fef451SLois Curfman McInnes 9034a2ae208SSatish Balay #undef __FUNCT__ 9044a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian" 905c2aafc4cSSatish Balay /*@C 906b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 907b4fd4287SBarry Smith provided context for evaluating the Jacobian. 908b4fd4287SBarry Smith 909c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 910c7afd0dbSLois Curfman McInnes 911b4fd4287SBarry Smith Input Parameter: 912b4fd4287SBarry Smith . snes - the nonlinear solver context 913b4fd4287SBarry Smith 914b4fd4287SBarry Smith Output Parameters: 915c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 916b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 91700036973SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 91800036973SBarry Smith - func - location to put Jacobian function (or PETSC_NULL) 919fee21e36SBarry Smith 92036851e7fSLois Curfman McInnes Level: advanced 92136851e7fSLois Curfman McInnes 922b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 923b4fd4287SBarry Smith @*/ 92463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 925b4fd4287SBarry Smith { 9263a40ed3dSBarry Smith PetscFunctionBegin; 9274482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 928b4fd4287SBarry Smith if (A) *A = snes->jacobian; 929b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 930b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 93100036973SBarry Smith if (func) *func = snes->computejacobian; 9323a40ed3dSBarry Smith PetscFunctionReturn(0); 933b4fd4287SBarry Smith } 934b4fd4287SBarry Smith 9359b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 93663dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 9379b94acceSBarry Smith 9384a2ae208SSatish Balay #undef __FUNCT__ 9394a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp" 9409b94acceSBarry Smith /*@ 9419b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 942272ac6f2SLois Curfman McInnes of a nonlinear solver. 9439b94acceSBarry Smith 944fee21e36SBarry Smith Collective on SNES 945fee21e36SBarry Smith 946c7afd0dbSLois Curfman McInnes Input Parameters: 947c7afd0dbSLois Curfman McInnes + snes - the SNES context 948c7afd0dbSLois Curfman McInnes - x - the solution vector 949c7afd0dbSLois Curfman McInnes 950272ac6f2SLois Curfman McInnes Notes: 951272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 952272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 953272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 954272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 955272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 956272ac6f2SLois Curfman McInnes 95736851e7fSLois Curfman McInnes Level: advanced 95836851e7fSLois Curfman McInnes 9599b94acceSBarry Smith .keywords: SNES, nonlinear, setup 9609b94acceSBarry Smith 9619b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 9629b94acceSBarry Smith @*/ 96363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes,Vec x) 9649b94acceSBarry Smith { 965dfbe8321SBarry Smith PetscErrorCode ierr; 9664b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 9673a40ed3dSBarry Smith 9683a40ed3dSBarry Smith PetscFunctionBegin; 9694482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 9704482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 971c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 9728ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 9739b94acceSBarry Smith 974b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 975c1f60f51SBarry Smith /* 976c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 977dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 978c1f60f51SBarry Smith */ 979c1f60f51SBarry Smith if (flg) { 980c1f60f51SBarry Smith Mat J; 9815a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 9825a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 98363ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator routines\n"));CHKERRQ(ierr); 9843a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 9853a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 986c1f60f51SBarry Smith } 98745fc7adcSBarry Smith 988c8d321e0SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) 98945fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 99045fc7adcSBarry Smith if (flg) { 99145fc7adcSBarry Smith Mat J; 99245fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 99345fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 99445fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 99545fc7adcSBarry Smith } 99632a4b47aSMatthew Knepley #endif 99745fc7adcSBarry Smith 998b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 999c1f60f51SBarry Smith /* 1000dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1001c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1002c1f60f51SBarry Smith */ 100331615d04SLois Curfman McInnes if (flg) { 1004272ac6f2SLois Curfman McInnes Mat J; 1005b5d62d44SBarry Smith KSP ksp; 100694b7f48cSBarry Smith PC pc; 1007f3ef73ceSBarry Smith 10085a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10093a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 101063ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"));CHKERRQ(ierr); 10113a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 10123a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 10133a7fca6bSBarry Smith 1014f3ef73ceSBarry Smith /* force no preconditioner */ 101594b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1016b5d62d44SBarry Smith ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1017a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1018a9815358SBarry Smith if (!flg) { 1019f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1020272ac6f2SLois Curfman McInnes } 1021a9815358SBarry Smith } 1022f3ef73ceSBarry Smith 102329bbc08cSBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 102429bbc08cSBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 102529bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1026a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 102729bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1028a8c6a408SBarry Smith } 1029a703fe33SLois Curfman McInnes 1030a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 10314b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 10326831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 103394b7f48cSBarry Smith KSP ksp; 103494b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1035186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1036a703fe33SLois Curfman McInnes } 10374b27c08aSLois Curfman McInnes 1038a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 103982bf6240SBarry Smith snes->setupcalled = 1; 10403a40ed3dSBarry Smith PetscFunctionReturn(0); 10419b94acceSBarry Smith } 10429b94acceSBarry Smith 10434a2ae208SSatish Balay #undef __FUNCT__ 10444a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 10459b94acceSBarry Smith /*@C 10469b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 10479b94acceSBarry Smith with SNESCreate(). 10489b94acceSBarry Smith 1049c7afd0dbSLois Curfman McInnes Collective on SNES 1050c7afd0dbSLois Curfman McInnes 10519b94acceSBarry Smith Input Parameter: 10529b94acceSBarry Smith . snes - the SNES context 10539b94acceSBarry Smith 105436851e7fSLois Curfman McInnes Level: beginner 105536851e7fSLois Curfman McInnes 10569b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 10579b94acceSBarry Smith 105863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 10599b94acceSBarry Smith @*/ 106063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 10619b94acceSBarry Smith { 106277431f27SBarry Smith PetscInt i; 10636849ba73SBarry Smith PetscErrorCode ierr; 10643a40ed3dSBarry Smith 10653a40ed3dSBarry Smith PetscFunctionBegin; 10664482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 10673a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1068d4bb536fSBarry Smith 1069be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 10700f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1071be0abb6dSBarry Smith 1072e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1073606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 10743a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 10753a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 10763ab0aad5SBarry Smith if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 107794b7f48cSBarry Smith ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1078522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1079b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1080b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1081b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1082b8a78c4aSBarry Smith } 1083b8a78c4aSBarry Smith } 1084a79aaaedSSatish Balay ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 10853a40ed3dSBarry Smith PetscFunctionReturn(0); 10869b94acceSBarry Smith } 10879b94acceSBarry Smith 10889b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 10899b94acceSBarry Smith 10904a2ae208SSatish Balay #undef __FUNCT__ 10914a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 10929b94acceSBarry Smith /*@ 1093d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 10949b94acceSBarry Smith 1095c7afd0dbSLois Curfman McInnes Collective on SNES 1096c7afd0dbSLois Curfman McInnes 10979b94acceSBarry Smith Input Parameters: 1098c7afd0dbSLois Curfman McInnes + snes - the SNES context 109970441072SBarry Smith . abstol - absolute convergence tolerance 110033174efeSLois Curfman McInnes . rtol - relative convergence tolerance 110133174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 110233174efeSLois Curfman McInnes of the change in the solution between steps 110333174efeSLois Curfman McInnes . maxit - maximum number of iterations 1104c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1105fee21e36SBarry Smith 110633174efeSLois Curfman McInnes Options Database Keys: 110770441072SBarry Smith + -snes_atol <abstol> - Sets abstol 1108c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1109c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1110c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1111c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 11129b94acceSBarry Smith 1113d7a720efSLois Curfman McInnes Notes: 11149b94acceSBarry Smith The default maximum number of iterations is 50. 11159b94acceSBarry Smith The default maximum number of function evaluations is 1000. 11169b94acceSBarry Smith 111736851e7fSLois Curfman McInnes Level: intermediate 111836851e7fSLois Curfman McInnes 111933174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 11209b94acceSBarry Smith 11212492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance() 11229b94acceSBarry Smith @*/ 112363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 11249b94acceSBarry Smith { 11253a40ed3dSBarry Smith PetscFunctionBegin; 11264482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 112770441072SBarry Smith if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1128d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1129d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1130d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1131d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11323a40ed3dSBarry Smith PetscFunctionReturn(0); 11339b94acceSBarry Smith } 11349b94acceSBarry Smith 11354a2ae208SSatish Balay #undef __FUNCT__ 11364a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11379b94acceSBarry Smith /*@ 113833174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 113933174efeSLois Curfman McInnes 1140c7afd0dbSLois Curfman McInnes Not Collective 1141c7afd0dbSLois Curfman McInnes 114233174efeSLois Curfman McInnes Input Parameters: 1143c7afd0dbSLois Curfman McInnes + snes - the SNES context 114470441072SBarry Smith . abstol - absolute convergence tolerance 114533174efeSLois Curfman McInnes . rtol - relative convergence tolerance 114633174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 114733174efeSLois Curfman McInnes of the change in the solution between steps 114833174efeSLois Curfman McInnes . maxit - maximum number of iterations 1149c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1150fee21e36SBarry Smith 115133174efeSLois Curfman McInnes Notes: 115233174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 115333174efeSLois Curfman McInnes 115436851e7fSLois Curfman McInnes Level: intermediate 115536851e7fSLois Curfman McInnes 115633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 115733174efeSLois Curfman McInnes 115833174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 115933174efeSLois Curfman McInnes @*/ 116063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 116133174efeSLois Curfman McInnes { 11623a40ed3dSBarry Smith PetscFunctionBegin; 11634482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 116470441072SBarry Smith if (abstol) *abstol = snes->abstol; 116533174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 116633174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 116733174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 116833174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 11693a40ed3dSBarry Smith PetscFunctionReturn(0); 117033174efeSLois Curfman McInnes } 117133174efeSLois Curfman McInnes 11724a2ae208SSatish Balay #undef __FUNCT__ 11734a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 117433174efeSLois Curfman McInnes /*@ 11759b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 11769b94acceSBarry Smith 1177fee21e36SBarry Smith Collective on SNES 1178fee21e36SBarry Smith 1179c7afd0dbSLois Curfman McInnes Input Parameters: 1180c7afd0dbSLois Curfman McInnes + snes - the SNES context 1181c7afd0dbSLois Curfman McInnes - tol - tolerance 1182c7afd0dbSLois Curfman McInnes 11839b94acceSBarry Smith Options Database Key: 1184c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 11859b94acceSBarry Smith 118636851e7fSLois Curfman McInnes Level: intermediate 118736851e7fSLois Curfman McInnes 11889b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 11899b94acceSBarry Smith 11902492ecdbSBarry Smith .seealso: SNESSetTolerances() 11919b94acceSBarry Smith @*/ 119263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 11939b94acceSBarry Smith { 11943a40ed3dSBarry Smith PetscFunctionBegin; 11954482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 11969b94acceSBarry Smith snes->deltatol = tol; 11973a40ed3dSBarry Smith PetscFunctionReturn(0); 11989b94acceSBarry Smith } 11999b94acceSBarry Smith 1200df9fa365SBarry Smith /* 1201df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1202df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1203df9fa365SBarry Smith macros instead of functions 1204df9fa365SBarry Smith */ 12054a2ae208SSatish Balay #undef __FUNCT__ 12064a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 120763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1208ce1608b8SBarry Smith { 1209dfbe8321SBarry Smith PetscErrorCode ierr; 1210ce1608b8SBarry Smith 1211ce1608b8SBarry Smith PetscFunctionBegin; 12124482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1213ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1214ce1608b8SBarry Smith PetscFunctionReturn(0); 1215ce1608b8SBarry Smith } 1216ce1608b8SBarry Smith 12174a2ae208SSatish Balay #undef __FUNCT__ 12184a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 121963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1220df9fa365SBarry Smith { 1221dfbe8321SBarry Smith PetscErrorCode ierr; 1222df9fa365SBarry Smith 1223df9fa365SBarry Smith PetscFunctionBegin; 1224df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1225df9fa365SBarry Smith PetscFunctionReturn(0); 1226df9fa365SBarry Smith } 1227df9fa365SBarry Smith 12284a2ae208SSatish Balay #undef __FUNCT__ 12294a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 123063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw) 1231df9fa365SBarry Smith { 1232dfbe8321SBarry Smith PetscErrorCode ierr; 1233df9fa365SBarry Smith 1234df9fa365SBarry Smith PetscFunctionBegin; 1235df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1236df9fa365SBarry Smith PetscFunctionReturn(0); 1237df9fa365SBarry Smith } 1238df9fa365SBarry Smith 12399b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12409b94acceSBarry Smith 12414a2ae208SSatish Balay #undef __FUNCT__ 12424a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12439b94acceSBarry Smith /*@C 12445cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12459b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12469b94acceSBarry Smith progress. 12479b94acceSBarry Smith 1248fee21e36SBarry Smith Collective on SNES 1249fee21e36SBarry Smith 1250c7afd0dbSLois Curfman McInnes Input Parameters: 1251c7afd0dbSLois Curfman McInnes + snes - the SNES context 1252c7afd0dbSLois Curfman McInnes . func - monitoring routine 1253b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1254b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desitre) 1255b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1256b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 12579b94acceSBarry Smith 1258c7afd0dbSLois Curfman McInnes Calling sequence of func: 1259a7cc72afSBarry Smith $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1260c7afd0dbSLois Curfman McInnes 1261c7afd0dbSLois Curfman McInnes + snes - the SNES context 1262c7afd0dbSLois Curfman McInnes . its - iteration number 1263c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 126440a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 12659b94acceSBarry Smith 12669665c990SLois Curfman McInnes Options Database Keys: 1267c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1268c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1269c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1270c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1271c7afd0dbSLois Curfman McInnes been hardwired into a code by 1272c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1273c7afd0dbSLois Curfman McInnes does not cancel those set via 1274c7afd0dbSLois Curfman McInnes the options database. 12759665c990SLois Curfman McInnes 1276639f9d9dSBarry Smith Notes: 12776bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 12786bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 12796bc08f3fSLois Curfman McInnes order in which they were set. 1280639f9d9dSBarry Smith 128136851e7fSLois Curfman McInnes Level: intermediate 128236851e7fSLois Curfman McInnes 12839b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 12849b94acceSBarry Smith 12855cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 12869b94acceSBarry Smith @*/ 128763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 12889b94acceSBarry Smith { 12893a40ed3dSBarry Smith PetscFunctionBegin; 12904482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1291639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 129229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1293639f9d9dSBarry Smith } 1294639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1295b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1296639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 12973a40ed3dSBarry Smith PetscFunctionReturn(0); 12989b94acceSBarry Smith } 12999b94acceSBarry Smith 13004a2ae208SSatish Balay #undef __FUNCT__ 13014a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 13025cd90555SBarry Smith /*@C 13035cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 13045cd90555SBarry Smith 1305c7afd0dbSLois Curfman McInnes Collective on SNES 1306c7afd0dbSLois Curfman McInnes 13075cd90555SBarry Smith Input Parameters: 13085cd90555SBarry Smith . snes - the SNES context 13095cd90555SBarry Smith 13101a480d89SAdministrator Options Database Key: 1311c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1312c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1313c7afd0dbSLois Curfman McInnes set via the options database 13145cd90555SBarry Smith 13155cd90555SBarry Smith Notes: 13165cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 13175cd90555SBarry Smith 131836851e7fSLois Curfman McInnes Level: intermediate 131936851e7fSLois Curfman McInnes 13205cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 13215cd90555SBarry Smith 13225cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 13235cd90555SBarry Smith @*/ 132463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes) 13255cd90555SBarry Smith { 13265cd90555SBarry Smith PetscFunctionBegin; 13274482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13285cd90555SBarry Smith snes->numbermonitors = 0; 13295cd90555SBarry Smith PetscFunctionReturn(0); 13305cd90555SBarry Smith } 13315cd90555SBarry Smith 13324a2ae208SSatish Balay #undef __FUNCT__ 13334a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13349b94acceSBarry Smith /*@C 13359b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13369b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13379b94acceSBarry Smith 1338fee21e36SBarry Smith Collective on SNES 1339fee21e36SBarry Smith 1340c7afd0dbSLois Curfman McInnes Input Parameters: 1341c7afd0dbSLois Curfman McInnes + snes - the SNES context 1342c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1343c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1344c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 13459b94acceSBarry Smith 1346c7afd0dbSLois Curfman McInnes Calling sequence of func: 13476849ba73SBarry Smith $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1348c7afd0dbSLois Curfman McInnes 1349c7afd0dbSLois Curfman McInnes + snes - the SNES context 1350c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1351184914b5SBarry Smith . reason - reason for convergence/divergence 1352c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 13534b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 13544b27c08aSLois Curfman McInnes - f - 2-norm of function 13559b94acceSBarry Smith 135636851e7fSLois Curfman McInnes Level: advanced 135736851e7fSLois Curfman McInnes 13589b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13599b94acceSBarry Smith 13604b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 13619b94acceSBarry Smith @*/ 136263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 13639b94acceSBarry Smith { 13643a40ed3dSBarry Smith PetscFunctionBegin; 13654482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13669b94acceSBarry Smith (snes)->converged = func; 13679b94acceSBarry Smith (snes)->cnvP = cctx; 13683a40ed3dSBarry Smith PetscFunctionReturn(0); 13699b94acceSBarry Smith } 13709b94acceSBarry Smith 13714a2ae208SSatish Balay #undef __FUNCT__ 13724a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 1373184914b5SBarry Smith /*@C 1374184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1375184914b5SBarry Smith 1376184914b5SBarry Smith Not Collective 1377184914b5SBarry Smith 1378184914b5SBarry Smith Input Parameter: 1379184914b5SBarry Smith . snes - the SNES context 1380184914b5SBarry Smith 1381184914b5SBarry Smith Output Parameter: 1382e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1383184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1384184914b5SBarry Smith 1385184914b5SBarry Smith Level: intermediate 1386184914b5SBarry Smith 1387184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1388184914b5SBarry Smith 1389184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1390184914b5SBarry Smith 13914b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1392184914b5SBarry Smith @*/ 139363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1394184914b5SBarry Smith { 1395184914b5SBarry Smith PetscFunctionBegin; 13964482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13974482741eSBarry Smith PetscValidPointer(reason,2); 1398184914b5SBarry Smith *reason = snes->reason; 1399184914b5SBarry Smith PetscFunctionReturn(0); 1400184914b5SBarry Smith } 1401184914b5SBarry Smith 14024a2ae208SSatish Balay #undef __FUNCT__ 14034a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1404c9005455SLois Curfman McInnes /*@ 1405c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1406c9005455SLois Curfman McInnes 1407fee21e36SBarry Smith Collective on SNES 1408fee21e36SBarry Smith 1409c7afd0dbSLois Curfman McInnes Input Parameters: 1410c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1411c7afd0dbSLois Curfman McInnes . a - array to hold history 1412cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1413758f92a0SBarry Smith . na - size of a and its 141464731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1415758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1416c7afd0dbSLois Curfman McInnes 1417c9005455SLois Curfman McInnes Notes: 14184b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1419c9005455SLois Curfman McInnes at each step. 1420c9005455SLois Curfman McInnes 1421c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1422c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1423c9005455SLois Curfman McInnes during the section of code that is being timed. 1424c9005455SLois Curfman McInnes 142536851e7fSLois Curfman McInnes Level: intermediate 142636851e7fSLois Curfman McInnes 1427c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1428758f92a0SBarry Smith 142908405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1430758f92a0SBarry Smith 1431c9005455SLois Curfman McInnes @*/ 143263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1433c9005455SLois Curfman McInnes { 14343a40ed3dSBarry Smith PetscFunctionBegin; 14354482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14364482741eSBarry Smith if (na) PetscValidScalarPointer(a,2); 1437c9005455SLois Curfman McInnes snes->conv_hist = a; 1438758f92a0SBarry Smith snes->conv_hist_its = its; 1439758f92a0SBarry Smith snes->conv_hist_max = na; 1440758f92a0SBarry Smith snes->conv_hist_reset = reset; 1441758f92a0SBarry Smith PetscFunctionReturn(0); 1442758f92a0SBarry Smith } 1443758f92a0SBarry Smith 14444a2ae208SSatish Balay #undef __FUNCT__ 14454a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 14460c4c9dddSBarry Smith /*@C 1447758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1448758f92a0SBarry Smith 1449758f92a0SBarry Smith Collective on SNES 1450758f92a0SBarry Smith 1451758f92a0SBarry Smith Input Parameter: 1452758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1453758f92a0SBarry Smith 1454758f92a0SBarry Smith Output Parameters: 1455758f92a0SBarry Smith . a - array to hold history 1456758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1457758f92a0SBarry Smith negative if not converged) for each solve. 1458758f92a0SBarry Smith - na - size of a and its 1459758f92a0SBarry Smith 1460758f92a0SBarry Smith Notes: 1461758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1462758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1463758f92a0SBarry Smith 1464758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1465758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1466758f92a0SBarry Smith during the section of code that is being timed. 1467758f92a0SBarry Smith 1468758f92a0SBarry Smith Level: intermediate 1469758f92a0SBarry Smith 1470758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1471758f92a0SBarry Smith 1472758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1473758f92a0SBarry Smith 1474758f92a0SBarry Smith @*/ 147563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1476758f92a0SBarry Smith { 1477758f92a0SBarry Smith PetscFunctionBegin; 14784482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1479758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1480758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1481758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 14823a40ed3dSBarry Smith PetscFunctionReturn(0); 1483c9005455SLois Curfman McInnes } 1484c9005455SLois Curfman McInnes 1485e74ef692SMatthew Knepley #undef __FUNCT__ 1486e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC" 1487ac226902SBarry Smith /*@C 148876b2cf59SMatthew Knepley SNESSetRhsBC - Sets the function which applies boundary conditions 148976b2cf59SMatthew Knepley to the Rhs of each system. 149076b2cf59SMatthew Knepley 149176b2cf59SMatthew Knepley Collective on SNES 149276b2cf59SMatthew Knepley 149376b2cf59SMatthew Knepley Input Parameters: 149476b2cf59SMatthew Knepley . snes - The nonlinear solver context 149576b2cf59SMatthew Knepley . func - The function 149676b2cf59SMatthew Knepley 149776b2cf59SMatthew Knepley Calling sequence of func: 149876b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx); 149976b2cf59SMatthew Knepley 150076b2cf59SMatthew Knepley . rhs - The current rhs vector 150176b2cf59SMatthew Knepley . ctx - The user-context 150276b2cf59SMatthew Knepley 150376b2cf59SMatthew Knepley Level: intermediate 150476b2cf59SMatthew Knepley 150576b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 150676b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate() 150776b2cf59SMatthew Knepley @*/ 150863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetRhsBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 150976b2cf59SMatthew Knepley { 151076b2cf59SMatthew Knepley PetscFunctionBegin; 15114482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 151276b2cf59SMatthew Knepley snes->applyrhsbc = func; 151376b2cf59SMatthew Knepley PetscFunctionReturn(0); 151476b2cf59SMatthew Knepley } 151576b2cf59SMatthew Knepley 1516e74ef692SMatthew Knepley #undef __FUNCT__ 1517e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC" 151876b2cf59SMatthew Knepley /*@ 151976b2cf59SMatthew Knepley SNESDefaultRhsBC - The default boundary condition function which does nothing. 152076b2cf59SMatthew Knepley 152176b2cf59SMatthew Knepley Not collective 152276b2cf59SMatthew Knepley 152376b2cf59SMatthew Knepley Input Parameters: 152476b2cf59SMatthew Knepley . snes - The nonlinear solver context 152576b2cf59SMatthew Knepley . rhs - The Rhs 152676b2cf59SMatthew Knepley . ctx - The user-context 152776b2cf59SMatthew Knepley 152876b2cf59SMatthew Knepley Level: beginner 152976b2cf59SMatthew Knepley 153076b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 153176b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate() 153276b2cf59SMatthew Knepley @*/ 153363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx) 153476b2cf59SMatthew Knepley { 153576b2cf59SMatthew Knepley PetscFunctionBegin; 153676b2cf59SMatthew Knepley PetscFunctionReturn(0); 153776b2cf59SMatthew Knepley } 153876b2cf59SMatthew Knepley 1539e74ef692SMatthew Knepley #undef __FUNCT__ 1540e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC" 1541ac226902SBarry Smith /*@C 154276b2cf59SMatthew Knepley SNESSetSolutionBC - Sets the function which applies boundary conditions 154376b2cf59SMatthew Knepley to the solution of each system. 154476b2cf59SMatthew Knepley 154576b2cf59SMatthew Knepley Collective on SNES 154676b2cf59SMatthew Knepley 154776b2cf59SMatthew Knepley Input Parameters: 154876b2cf59SMatthew Knepley . snes - The nonlinear solver context 154976b2cf59SMatthew Knepley . func - The function 155076b2cf59SMatthew Knepley 155176b2cf59SMatthew Knepley Calling sequence of func: 15529d06fb6bSMatthew Knepley . func (SNES snes, Vec rsol, void *ctx); 155376b2cf59SMatthew Knepley 155476b2cf59SMatthew Knepley . sol - The current solution vector 155576b2cf59SMatthew Knepley . ctx - The user-context 155676b2cf59SMatthew Knepley 155776b2cf59SMatthew Knepley Level: intermediate 155876b2cf59SMatthew Knepley 155976b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 156076b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate() 156176b2cf59SMatthew Knepley @*/ 156263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolutionBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 156376b2cf59SMatthew Knepley { 156476b2cf59SMatthew Knepley PetscFunctionBegin; 15654482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 156676b2cf59SMatthew Knepley snes->applysolbc = func; 156776b2cf59SMatthew Knepley PetscFunctionReturn(0); 156876b2cf59SMatthew Knepley } 156976b2cf59SMatthew Knepley 1570e74ef692SMatthew Knepley #undef __FUNCT__ 1571e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC" 157276b2cf59SMatthew Knepley /*@ 157376b2cf59SMatthew Knepley SNESDefaultSolutionBC - The default boundary condition function which does nothing. 157476b2cf59SMatthew Knepley 157576b2cf59SMatthew Knepley Not collective 157676b2cf59SMatthew Knepley 157776b2cf59SMatthew Knepley Input Parameters: 157876b2cf59SMatthew Knepley . snes - The nonlinear solver context 157976b2cf59SMatthew Knepley . sol - The solution 158076b2cf59SMatthew Knepley . ctx - The user-context 158176b2cf59SMatthew Knepley 158276b2cf59SMatthew Knepley Level: beginner 158376b2cf59SMatthew Knepley 158476b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 158576b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate() 158676b2cf59SMatthew Knepley @*/ 158763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx) 158876b2cf59SMatthew Knepley { 158976b2cf59SMatthew Knepley PetscFunctionBegin; 159076b2cf59SMatthew Knepley PetscFunctionReturn(0); 159176b2cf59SMatthew Knepley } 159276b2cf59SMatthew Knepley 1593e74ef692SMatthew Knepley #undef __FUNCT__ 1594e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 1595ac226902SBarry Smith /*@C 159676b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 159776b2cf59SMatthew Knepley at the beginning of every step of the iteration. 159876b2cf59SMatthew Knepley 159976b2cf59SMatthew Knepley Collective on SNES 160076b2cf59SMatthew Knepley 160176b2cf59SMatthew Knepley Input Parameters: 160276b2cf59SMatthew Knepley . snes - The nonlinear solver context 160376b2cf59SMatthew Knepley . func - The function 160476b2cf59SMatthew Knepley 160576b2cf59SMatthew Knepley Calling sequence of func: 1606b5d30489SBarry Smith . func (SNES snes, PetscInt step); 160776b2cf59SMatthew Knepley 160876b2cf59SMatthew Knepley . step - The current step of the iteration 160976b2cf59SMatthew Knepley 161076b2cf59SMatthew Knepley Level: intermediate 161176b2cf59SMatthew Knepley 161276b2cf59SMatthew Knepley .keywords: SNES, update 1613b5d30489SBarry Smith 161476b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 161576b2cf59SMatthew Knepley @*/ 161663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 161776b2cf59SMatthew Knepley { 161876b2cf59SMatthew Knepley PetscFunctionBegin; 16194482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 162076b2cf59SMatthew Knepley snes->update = func; 162176b2cf59SMatthew Knepley PetscFunctionReturn(0); 162276b2cf59SMatthew Knepley } 162376b2cf59SMatthew Knepley 1624e74ef692SMatthew Knepley #undef __FUNCT__ 1625e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 162676b2cf59SMatthew Knepley /*@ 162776b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 162876b2cf59SMatthew Knepley 162976b2cf59SMatthew Knepley Not collective 163076b2cf59SMatthew Knepley 163176b2cf59SMatthew Knepley Input Parameters: 163276b2cf59SMatthew Knepley . snes - The nonlinear solver context 163376b2cf59SMatthew Knepley . step - The current step of the iteration 163476b2cf59SMatthew Knepley 1635205452f4SMatthew Knepley Level: intermediate 1636205452f4SMatthew Knepley 163776b2cf59SMatthew Knepley .keywords: SNES, update 163876b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 163976b2cf59SMatthew Knepley @*/ 164063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 164176b2cf59SMatthew Knepley { 164276b2cf59SMatthew Knepley PetscFunctionBegin; 164376b2cf59SMatthew Knepley PetscFunctionReturn(0); 164476b2cf59SMatthew Knepley } 164576b2cf59SMatthew Knepley 16464a2ae208SSatish Balay #undef __FUNCT__ 16474a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 16489b94acceSBarry Smith /* 16499b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 16509b94acceSBarry Smith positive parameter delta. 16519b94acceSBarry Smith 16529b94acceSBarry Smith Input Parameters: 1653c7afd0dbSLois Curfman McInnes + snes - the SNES context 16549b94acceSBarry Smith . y - approximate solution of linear system 16559b94acceSBarry Smith . fnorm - 2-norm of current function 1656c7afd0dbSLois Curfman McInnes - delta - trust region size 16579b94acceSBarry Smith 16589b94acceSBarry Smith Output Parameters: 1659c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 16609b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 16619b94acceSBarry Smith region, and exceeds zero otherwise. 1662c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 16639b94acceSBarry Smith 16649b94acceSBarry Smith Note: 16654b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 16669b94acceSBarry Smith is set to be the maximum allowable step size. 16679b94acceSBarry Smith 16689b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16699b94acceSBarry Smith */ 1670dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 16719b94acceSBarry Smith { 1672064f8208SBarry Smith PetscReal nrm; 1673ea709b57SSatish Balay PetscScalar cnorm; 1674dfbe8321SBarry Smith PetscErrorCode ierr; 16753a40ed3dSBarry Smith 16763a40ed3dSBarry Smith PetscFunctionBegin; 16774482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16784482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1679c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,2); 1680184914b5SBarry Smith 1681064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1682064f8208SBarry Smith if (nrm > *delta) { 1683064f8208SBarry Smith nrm = *delta/nrm; 1684064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1685064f8208SBarry Smith cnorm = nrm; 1686329f5518SBarry Smith ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 16879b94acceSBarry Smith *ynorm = *delta; 16889b94acceSBarry Smith } else { 16899b94acceSBarry Smith *gpnorm = 0.0; 1690064f8208SBarry Smith *ynorm = nrm; 16919b94acceSBarry Smith } 16923a40ed3dSBarry Smith PetscFunctionReturn(0); 16939b94acceSBarry Smith } 16949b94acceSBarry Smith 16955968eb51SBarry Smith static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero", 16965968eb51SBarry Smith "not currently used", 16975968eb51SBarry Smith "line search failed", 16985968eb51SBarry Smith "reach maximum number of iterations", 16995968eb51SBarry Smith "function norm became NaN (not a number)", 17005968eb51SBarry Smith "not currently used", 17015968eb51SBarry Smith "number of function computations exceeded", 1702*d5e45103SBarry Smith "arguments to user function are not in domain", 17035968eb51SBarry Smith "still iterating", 17045968eb51SBarry Smith "not currently used", 17055968eb51SBarry Smith "absolute size of function norm", 17065968eb51SBarry Smith "relative decrease in function norm", 17075968eb51SBarry Smith "step size is small", 17085968eb51SBarry Smith "not currently used", 17095968eb51SBarry Smith "not currently used", 17105968eb51SBarry Smith "small trust region"}; 17115968eb51SBarry Smith 17124a2ae208SSatish Balay #undef __FUNCT__ 17134a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 17149b94acceSBarry Smith /*@ 17159b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 171663a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 17179b94acceSBarry Smith 1718c7afd0dbSLois Curfman McInnes Collective on SNES 1719c7afd0dbSLois Curfman McInnes 1720b2002411SLois Curfman McInnes Input Parameters: 1721c7afd0dbSLois Curfman McInnes + snes - the SNES context 1722c7afd0dbSLois Curfman McInnes - x - the solution vector 17239b94acceSBarry Smith 1724b2002411SLois Curfman McInnes Notes: 17258ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 17268ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 17278ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 17288ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 17298ddd3da0SLois Curfman McInnes 173036851e7fSLois Curfman McInnes Level: beginner 173136851e7fSLois Curfman McInnes 17329b94acceSBarry Smith .keywords: SNES, nonlinear, solve 17339b94acceSBarry Smith 17343ab0aad5SBarry Smith .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs() 17359b94acceSBarry Smith @*/ 173663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec x) 17379b94acceSBarry Smith { 1738dfbe8321SBarry Smith PetscErrorCode ierr; 1739f1af5d2fSBarry Smith PetscTruth flg; 1740052efed2SBarry Smith 17413a40ed3dSBarry Smith PetscFunctionBegin; 17424482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17434482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1744c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 17451302d50aSBarry Smith if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 174674637425SBarry Smith 174782bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1748c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 1749abc0a331SBarry Smith if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1750d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 175150ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1752*d5e45103SBarry Smith 1753*d5e45103SBarry Smith ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN); 1754*d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 1755*d5e45103SBarry Smith /* this means that a caller above me has also tryed this exception so I don't handle it here */ 175619717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr); 1757*d5e45103SBarry Smith } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) { 1758*d5e45103SBarry Smith /* translate exception into SNES not converged reason */ 1759*d5e45103SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1760*d5e45103SBarry Smith } else CHKERRQ(ierr); 1761*d5e45103SBarry Smith 1762d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1763b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 17648b179ff8SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1765da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1766da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 17675968eb51SBarry Smith if (snes->printreason) { 17685968eb51SBarry Smith if (snes->reason > 0) { 17695968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 17705968eb51SBarry Smith } else { 17715968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 17725968eb51SBarry Smith } 17735968eb51SBarry Smith } 17745968eb51SBarry Smith 17753a40ed3dSBarry Smith PetscFunctionReturn(0); 17769b94acceSBarry Smith } 17779b94acceSBarry Smith 17789b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 17799b94acceSBarry Smith 17804a2ae208SSatish Balay #undef __FUNCT__ 17814a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 178282bf6240SBarry Smith /*@C 17834b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 17849b94acceSBarry Smith 1785fee21e36SBarry Smith Collective on SNES 1786fee21e36SBarry Smith 1787c7afd0dbSLois Curfman McInnes Input Parameters: 1788c7afd0dbSLois Curfman McInnes + snes - the SNES context 1789454a90a3SBarry Smith - type - a known method 1790c7afd0dbSLois Curfman McInnes 1791c7afd0dbSLois Curfman McInnes Options Database Key: 1792454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1793c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1794ae12b187SLois Curfman McInnes 17959b94acceSBarry Smith Notes: 1796e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 17974b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1798c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17994b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1800c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 18019b94acceSBarry Smith 1802ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1803ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1804ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1805ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1806ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1807ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1808ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1809ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1810ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1811b0a32e0cSBarry Smith appropriate method. 181236851e7fSLois Curfman McInnes 181336851e7fSLois Curfman McInnes Level: intermediate 1814a703fe33SLois Curfman McInnes 1815454a90a3SBarry Smith .keywords: SNES, set, type 1816435da068SBarry Smith 1817435da068SBarry Smith .seealso: SNESType, SNESCreate() 1818435da068SBarry Smith 18199b94acceSBarry Smith @*/ 182063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type) 18219b94acceSBarry Smith { 1822dfbe8321SBarry Smith PetscErrorCode ierr,(*r)(SNES); 18236831982aSBarry Smith PetscTruth match; 18243a40ed3dSBarry Smith 18253a40ed3dSBarry Smith PetscFunctionBegin; 18264482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18274482741eSBarry Smith PetscValidCharPointer(type,2); 182882bf6240SBarry Smith 18296831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 18300f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 183182bf6240SBarry Smith 183282bf6240SBarry Smith if (snes->setupcalled) { 1833e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 183482bf6240SBarry Smith snes->data = 0; 183582bf6240SBarry Smith } 18367d1a2b2bSBarry Smith 18379b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 183882bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1839b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1840958c9bccSBarry Smith if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1841606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 184282bf6240SBarry Smith snes->data = 0; 18433a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1844454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 18453a40ed3dSBarry Smith PetscFunctionReturn(0); 18469b94acceSBarry Smith } 18479b94acceSBarry Smith 1848a847f771SSatish Balay 18499b94acceSBarry Smith /* --------------------------------------------------------------------- */ 18504a2ae208SSatish Balay #undef __FUNCT__ 18514a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 18529b94acceSBarry Smith /*@C 18539b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1854f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 18559b94acceSBarry Smith 1856fee21e36SBarry Smith Not Collective 1857fee21e36SBarry Smith 185836851e7fSLois Curfman McInnes Level: advanced 185936851e7fSLois Curfman McInnes 18609b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 18619b94acceSBarry Smith 18629b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 18639b94acceSBarry Smith @*/ 186463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 18659b94acceSBarry Smith { 1866dfbe8321SBarry Smith PetscErrorCode ierr; 186782bf6240SBarry Smith 18683a40ed3dSBarry Smith PetscFunctionBegin; 186982bf6240SBarry Smith if (SNESList) { 1870b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 187182bf6240SBarry Smith SNESList = 0; 18729b94acceSBarry Smith } 18734c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 18743a40ed3dSBarry Smith PetscFunctionReturn(0); 18759b94acceSBarry Smith } 18769b94acceSBarry Smith 18774a2ae208SSatish Balay #undef __FUNCT__ 18784a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 18799b94acceSBarry Smith /*@C 18809a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 18819b94acceSBarry Smith 1882c7afd0dbSLois Curfman McInnes Not Collective 1883c7afd0dbSLois Curfman McInnes 18849b94acceSBarry Smith Input Parameter: 18854b0e389bSBarry Smith . snes - nonlinear solver context 18869b94acceSBarry Smith 18879b94acceSBarry Smith Output Parameter: 18883a7fca6bSBarry Smith . type - SNES method (a character string) 18899b94acceSBarry Smith 189036851e7fSLois Curfman McInnes Level: intermediate 189136851e7fSLois Curfman McInnes 1892454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 18939b94acceSBarry Smith @*/ 189463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type) 18959b94acceSBarry Smith { 18963a40ed3dSBarry Smith PetscFunctionBegin; 18974482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18984482741eSBarry Smith PetscValidPointer(type,2); 1899454a90a3SBarry Smith *type = snes->type_name; 19003a40ed3dSBarry Smith PetscFunctionReturn(0); 19019b94acceSBarry Smith } 19029b94acceSBarry Smith 19034a2ae208SSatish Balay #undef __FUNCT__ 19044a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 19059b94acceSBarry Smith /*@C 19069b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 19079b94acceSBarry Smith stored. 19089b94acceSBarry Smith 1909c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1910c7afd0dbSLois Curfman McInnes 19119b94acceSBarry Smith Input Parameter: 19129b94acceSBarry Smith . snes - the SNES context 19139b94acceSBarry Smith 19149b94acceSBarry Smith Output Parameter: 19159b94acceSBarry Smith . x - the solution 19169b94acceSBarry Smith 191736851e7fSLois Curfman McInnes Level: advanced 191836851e7fSLois Curfman McInnes 19199b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 19209b94acceSBarry Smith 19214b27c08aSLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetSolutionUpdate() 19229b94acceSBarry Smith @*/ 192363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 19249b94acceSBarry Smith { 19253a40ed3dSBarry Smith PetscFunctionBegin; 19264482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19274482741eSBarry Smith PetscValidPointer(x,2); 19289b94acceSBarry Smith *x = snes->vec_sol_always; 19293a40ed3dSBarry Smith PetscFunctionReturn(0); 19309b94acceSBarry Smith } 19319b94acceSBarry Smith 19324a2ae208SSatish Balay #undef __FUNCT__ 19334a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 19349b94acceSBarry Smith /*@C 19359b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 19369b94acceSBarry Smith stored. 19379b94acceSBarry Smith 1938c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1939c7afd0dbSLois Curfman McInnes 19409b94acceSBarry Smith Input Parameter: 19419b94acceSBarry Smith . snes - the SNES context 19429b94acceSBarry Smith 19439b94acceSBarry Smith Output Parameter: 19449b94acceSBarry Smith . x - the solution update 19459b94acceSBarry Smith 194636851e7fSLois Curfman McInnes Level: advanced 194736851e7fSLois Curfman McInnes 19489b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 19499b94acceSBarry Smith 19509b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 19519b94acceSBarry Smith @*/ 195263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 19539b94acceSBarry Smith { 19543a40ed3dSBarry Smith PetscFunctionBegin; 19554482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19564482741eSBarry Smith PetscValidPointer(x,2); 19579b94acceSBarry Smith *x = snes->vec_sol_update_always; 19583a40ed3dSBarry Smith PetscFunctionReturn(0); 19599b94acceSBarry Smith } 19609b94acceSBarry Smith 19614a2ae208SSatish Balay #undef __FUNCT__ 19624a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 19639b94acceSBarry Smith /*@C 19643638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 19659b94acceSBarry Smith 1966c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1967c7afd0dbSLois Curfman McInnes 19689b94acceSBarry Smith Input Parameter: 19699b94acceSBarry Smith . snes - the SNES context 19709b94acceSBarry Smith 19719b94acceSBarry Smith Output Parameter: 19727bf4e008SBarry Smith + r - the function (or PETSC_NULL) 197300036973SBarry Smith . ctx - the function context (or PETSC_NULL) 197400036973SBarry Smith - func - the function (or PETSC_NULL) 19759b94acceSBarry Smith 197636851e7fSLois Curfman McInnes Level: advanced 197736851e7fSLois Curfman McInnes 1978a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 19799b94acceSBarry Smith 19804b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 19819b94acceSBarry Smith @*/ 198263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,void **ctx,PetscErrorCode (**func)(SNES,Vec,Vec,void*)) 19839b94acceSBarry Smith { 19843a40ed3dSBarry Smith PetscFunctionBegin; 19854482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19867bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 19877bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 198800036973SBarry Smith if (func) *func = snes->computefunction; 19893a40ed3dSBarry Smith PetscFunctionReturn(0); 19909b94acceSBarry Smith } 19919b94acceSBarry Smith 19924a2ae208SSatish Balay #undef __FUNCT__ 19934a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 19943c7409f5SSatish Balay /*@C 19953c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1996d850072dSLois Curfman McInnes SNES options in the database. 19973c7409f5SSatish Balay 1998fee21e36SBarry Smith Collective on SNES 1999fee21e36SBarry Smith 2000c7afd0dbSLois Curfman McInnes Input Parameter: 2001c7afd0dbSLois Curfman McInnes + snes - the SNES context 2002c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2003c7afd0dbSLois Curfman McInnes 2004d850072dSLois Curfman McInnes Notes: 2005a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2006c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2007d850072dSLois Curfman McInnes 200836851e7fSLois Curfman McInnes Level: advanced 200936851e7fSLois Curfman McInnes 20103c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 2011a86d99e1SLois Curfman McInnes 2012a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 20133c7409f5SSatish Balay @*/ 201463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 20153c7409f5SSatish Balay { 2016dfbe8321SBarry Smith PetscErrorCode ierr; 20173c7409f5SSatish Balay 20183a40ed3dSBarry Smith PetscFunctionBegin; 20194482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2020639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 202194b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20223a40ed3dSBarry Smith PetscFunctionReturn(0); 20233c7409f5SSatish Balay } 20243c7409f5SSatish Balay 20254a2ae208SSatish Balay #undef __FUNCT__ 20264a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 20273c7409f5SSatish Balay /*@C 2028f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2029d850072dSLois Curfman McInnes SNES options in the database. 20303c7409f5SSatish Balay 2031fee21e36SBarry Smith Collective on SNES 2032fee21e36SBarry Smith 2033c7afd0dbSLois Curfman McInnes Input Parameters: 2034c7afd0dbSLois Curfman McInnes + snes - the SNES context 2035c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2036c7afd0dbSLois Curfman McInnes 2037d850072dSLois Curfman McInnes Notes: 2038a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2039c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2040d850072dSLois Curfman McInnes 204136851e7fSLois Curfman McInnes Level: advanced 204236851e7fSLois Curfman McInnes 20433c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2044a86d99e1SLois Curfman McInnes 2045a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20463c7409f5SSatish Balay @*/ 204763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 20483c7409f5SSatish Balay { 2049dfbe8321SBarry Smith PetscErrorCode ierr; 20503c7409f5SSatish Balay 20513a40ed3dSBarry Smith PetscFunctionBegin; 20524482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2053639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 205494b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20553a40ed3dSBarry Smith PetscFunctionReturn(0); 20563c7409f5SSatish Balay } 20573c7409f5SSatish Balay 20584a2ae208SSatish Balay #undef __FUNCT__ 20594a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 20609ab63eb5SSatish Balay /*@C 20613c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20623c7409f5SSatish Balay SNES options in the database. 20633c7409f5SSatish Balay 2064c7afd0dbSLois Curfman McInnes Not Collective 2065c7afd0dbSLois Curfman McInnes 20663c7409f5SSatish Balay Input Parameter: 20673c7409f5SSatish Balay . snes - the SNES context 20683c7409f5SSatish Balay 20693c7409f5SSatish Balay Output Parameter: 20703c7409f5SSatish Balay . prefix - pointer to the prefix string used 20713c7409f5SSatish Balay 20729ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 20739ab63eb5SSatish Balay sufficient length to hold the prefix. 20749ab63eb5SSatish Balay 207536851e7fSLois Curfman McInnes Level: advanced 207636851e7fSLois Curfman McInnes 20773c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2078a86d99e1SLois Curfman McInnes 2079a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20803c7409f5SSatish Balay @*/ 208163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,char *prefix[]) 20823c7409f5SSatish Balay { 2083dfbe8321SBarry Smith PetscErrorCode ierr; 20843c7409f5SSatish Balay 20853a40ed3dSBarry Smith PetscFunctionBegin; 20864482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2087639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 20883a40ed3dSBarry Smith PetscFunctionReturn(0); 20893c7409f5SSatish Balay } 20903c7409f5SSatish Balay 2091b2002411SLois Curfman McInnes 20924a2ae208SSatish Balay #undef __FUNCT__ 20934a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 20943cea93caSBarry Smith /*@C 20953cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 20963cea93caSBarry Smith 20977f6c08e0SMatthew Knepley Level: advanced 20983cea93caSBarry Smith @*/ 209963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2100b2002411SLois Curfman McInnes { 2101e2d1d2b7SBarry Smith char fullname[PETSC_MAX_PATH_LEN]; 2102dfbe8321SBarry Smith PetscErrorCode ierr; 2103b2002411SLois Curfman McInnes 2104b2002411SLois Curfman McInnes PetscFunctionBegin; 2105b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2106c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2107b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2108b2002411SLois Curfman McInnes } 2109da9b6338SBarry Smith 2110da9b6338SBarry Smith #undef __FUNCT__ 2111da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 211263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2113da9b6338SBarry Smith { 2114dfbe8321SBarry Smith PetscErrorCode ierr; 211577431f27SBarry Smith PetscInt N,i,j; 2116da9b6338SBarry Smith Vec u,uh,fh; 2117da9b6338SBarry Smith PetscScalar value; 2118da9b6338SBarry Smith PetscReal norm; 2119da9b6338SBarry Smith 2120da9b6338SBarry Smith PetscFunctionBegin; 2121da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2122da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2123da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2124da9b6338SBarry Smith 2125da9b6338SBarry Smith /* currently only works for sequential */ 2126da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2127da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2128da9b6338SBarry Smith for (i=0; i<N; i++) { 2129da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 213077431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2131da9b6338SBarry Smith for (j=-10; j<11; j++) { 2132ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2133da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 21343ab0aad5SBarry Smith ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2135da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 213677431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2137da9b6338SBarry Smith value = -value; 2138da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2139da9b6338SBarry Smith } 2140da9b6338SBarry Smith } 2141da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2142da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2143da9b6338SBarry Smith PetscFunctionReturn(0); 2144da9b6338SBarry Smith } 2145