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; 588*ed1caa07SMatthew Knepley PetscValidPointer(outsnes,2); 5898ba1e511SMatthew Knepley *outsnes = PETSC_NULL; 5908ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 5918ba1e511SMatthew Knepley ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 5928ba1e511SMatthew Knepley #endif 5938ba1e511SMatthew Knepley 59452e6d16bSBarry Smith ierr = PetscHeaderCreate(snes,_p_SNES,PetscInt,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 595e24b481bSBarry Smith snes->bops->publish = SNESPublish_Petsc; 5969b94acceSBarry Smith snes->max_its = 50; 5979750a799SBarry Smith snes->max_funcs = 10000; 5989b94acceSBarry Smith snes->norm = 0.0; 599b4874afaSBarry Smith snes->rtol = 1.e-8; 600b4874afaSBarry Smith snes->ttol = 0.0; 60170441072SBarry Smith snes->abstol = 1.e-50; 6029b94acceSBarry Smith snes->xtol = 1.e-8; 6034b27c08aSLois Curfman McInnes snes->deltatol = 1.e-12; 6049b94acceSBarry Smith snes->nfuncs = 0; 60550ffb88aSMatthew Knepley snes->numFailures = 0; 60650ffb88aSMatthew Knepley snes->maxFailures = 1; 6077a00f4a9SLois Curfman McInnes snes->linear_its = 0; 608639f9d9dSBarry Smith snes->numbermonitors = 0; 6099b94acceSBarry Smith snes->data = 0; 6109b94acceSBarry Smith snes->view = 0; 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; 770d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 77119717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 77219717074SBarry Smith } 773d5e45103SBarry 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) 91770e92668SMatthew Knepley . func - location to put Jacobian function (or PETSC_NULL) 91870e92668SMatthew Knepley - ctx - location to stash Jacobian ctx (or PETSC_NULL) 919fee21e36SBarry Smith 92036851e7fSLois Curfman McInnes Level: advanced 92136851e7fSLois Curfman McInnes 922b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 923b4fd4287SBarry Smith @*/ 92470e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 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; 93000036973SBarry Smith if (func) *func = snes->computejacobian; 93170e92668SMatthew Knepley if (ctx) *ctx = snes->jacP; 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: 94770e92668SMatthew Knepley . snes - the SNES context 948c7afd0dbSLois Curfman McInnes 949272ac6f2SLois Curfman McInnes Notes: 950272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 951272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 952272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 953272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 954272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 955272ac6f2SLois Curfman McInnes 95636851e7fSLois Curfman McInnes Level: advanced 95736851e7fSLois Curfman McInnes 9589b94acceSBarry Smith .keywords: SNES, nonlinear, setup 9599b94acceSBarry Smith 9609b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 9619b94acceSBarry Smith @*/ 96270e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 9639b94acceSBarry Smith { 964dfbe8321SBarry Smith PetscErrorCode ierr; 9654b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 9663a40ed3dSBarry Smith 9673a40ed3dSBarry Smith PetscFunctionBegin; 9684482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 9699b94acceSBarry Smith 970b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 971c1f60f51SBarry Smith /* 972c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 973dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 974c1f60f51SBarry Smith */ 975c1f60f51SBarry Smith if (flg) { 976c1f60f51SBarry Smith Mat J; 9775a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 9785a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 97963ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator routines\n"));CHKERRQ(ierr); 9803a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 9813a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 982c1f60f51SBarry Smith } 98345fc7adcSBarry Smith 984c8d321e0SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) 98545fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 98645fc7adcSBarry Smith if (flg) { 98745fc7adcSBarry Smith Mat J; 98845fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 98945fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 99045fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 99145fc7adcSBarry Smith } 99232a4b47aSMatthew Knepley #endif 99345fc7adcSBarry Smith 994b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 995c1f60f51SBarry Smith /* 996dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 997c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 998c1f60f51SBarry Smith */ 99931615d04SLois Curfman McInnes if (flg) { 1000272ac6f2SLois Curfman McInnes Mat J; 1001b5d62d44SBarry Smith KSP ksp; 100294b7f48cSBarry Smith PC pc; 1003f3ef73ceSBarry Smith 10045a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10053a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 100663ba0a88SBarry Smith ierr = PetscLogInfo((snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"));CHKERRQ(ierr); 10073a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 10083a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 10093a7fca6bSBarry Smith 1010f3ef73ceSBarry Smith /* force no preconditioner */ 101194b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1012b5d62d44SBarry Smith ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1013a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1014a9815358SBarry Smith if (!flg) { 1015f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1016272ac6f2SLois Curfman McInnes } 1017a9815358SBarry Smith } 1018f3ef73ceSBarry Smith 101929bbc08cSBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 102029bbc08cSBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 102129bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1022a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 102329bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1024a8c6a408SBarry Smith } 1025a703fe33SLois Curfman McInnes 1026a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 10274b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 10286831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 102994b7f48cSBarry Smith KSP ksp; 103094b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1031186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1032a703fe33SLois Curfman McInnes } 10334b27c08aSLois Curfman McInnes 1034a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 103582bf6240SBarry Smith snes->setupcalled = 1; 10363a40ed3dSBarry Smith PetscFunctionReturn(0); 10379b94acceSBarry Smith } 10389b94acceSBarry Smith 10394a2ae208SSatish Balay #undef __FUNCT__ 10404a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 10419b94acceSBarry Smith /*@C 10429b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 10439b94acceSBarry Smith with SNESCreate(). 10449b94acceSBarry Smith 1045c7afd0dbSLois Curfman McInnes Collective on SNES 1046c7afd0dbSLois Curfman McInnes 10479b94acceSBarry Smith Input Parameter: 10489b94acceSBarry Smith . snes - the SNES context 10499b94acceSBarry Smith 105036851e7fSLois Curfman McInnes Level: beginner 105136851e7fSLois Curfman McInnes 10529b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 10539b94acceSBarry Smith 105463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 10559b94acceSBarry Smith @*/ 105663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 10579b94acceSBarry Smith { 105877431f27SBarry Smith PetscInt i; 10596849ba73SBarry Smith PetscErrorCode ierr; 10603a40ed3dSBarry Smith 10613a40ed3dSBarry Smith PetscFunctionBegin; 10624482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 10633a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1064d4bb536fSBarry Smith 1065be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 10660f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1067be0abb6dSBarry Smith 1068e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1069606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 10703a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 10713a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 10723ab0aad5SBarry Smith if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 107394b7f48cSBarry Smith ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1074522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1075b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1076b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1077b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1078b8a78c4aSBarry Smith } 1079b8a78c4aSBarry Smith } 1080a79aaaedSSatish Balay ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 10813a40ed3dSBarry Smith PetscFunctionReturn(0); 10829b94acceSBarry Smith } 10839b94acceSBarry Smith 10849b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 10859b94acceSBarry Smith 10864a2ae208SSatish Balay #undef __FUNCT__ 10874a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 10889b94acceSBarry Smith /*@ 1089d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 10909b94acceSBarry Smith 1091c7afd0dbSLois Curfman McInnes Collective on SNES 1092c7afd0dbSLois Curfman McInnes 10939b94acceSBarry Smith Input Parameters: 1094c7afd0dbSLois Curfman McInnes + snes - the SNES context 109570441072SBarry Smith . abstol - absolute convergence tolerance 109633174efeSLois Curfman McInnes . rtol - relative convergence tolerance 109733174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 109833174efeSLois Curfman McInnes of the change in the solution between steps 109933174efeSLois Curfman McInnes . maxit - maximum number of iterations 1100c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1101fee21e36SBarry Smith 110233174efeSLois Curfman McInnes Options Database Keys: 110370441072SBarry Smith + -snes_atol <abstol> - Sets abstol 1104c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1105c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1106c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1107c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 11089b94acceSBarry Smith 1109d7a720efSLois Curfman McInnes Notes: 11109b94acceSBarry Smith The default maximum number of iterations is 50. 11119b94acceSBarry Smith The default maximum number of function evaluations is 1000. 11129b94acceSBarry Smith 111336851e7fSLois Curfman McInnes Level: intermediate 111436851e7fSLois Curfman McInnes 111533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 11169b94acceSBarry Smith 11172492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance() 11189b94acceSBarry Smith @*/ 111963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 11209b94acceSBarry Smith { 11213a40ed3dSBarry Smith PetscFunctionBegin; 11224482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 112370441072SBarry Smith if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1124d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1125d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1126d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1127d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11283a40ed3dSBarry Smith PetscFunctionReturn(0); 11299b94acceSBarry Smith } 11309b94acceSBarry Smith 11314a2ae208SSatish Balay #undef __FUNCT__ 11324a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11339b94acceSBarry Smith /*@ 113433174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 113533174efeSLois Curfman McInnes 1136c7afd0dbSLois Curfman McInnes Not Collective 1137c7afd0dbSLois Curfman McInnes 113833174efeSLois Curfman McInnes Input Parameters: 1139c7afd0dbSLois Curfman McInnes + snes - the SNES context 114070441072SBarry Smith . abstol - absolute convergence tolerance 114133174efeSLois Curfman McInnes . rtol - relative convergence tolerance 114233174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 114333174efeSLois Curfman McInnes of the change in the solution between steps 114433174efeSLois Curfman McInnes . maxit - maximum number of iterations 1145c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1146fee21e36SBarry Smith 114733174efeSLois Curfman McInnes Notes: 114833174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 114933174efeSLois Curfman McInnes 115036851e7fSLois Curfman McInnes Level: intermediate 115136851e7fSLois Curfman McInnes 115233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 115333174efeSLois Curfman McInnes 115433174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 115533174efeSLois Curfman McInnes @*/ 115663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 115733174efeSLois Curfman McInnes { 11583a40ed3dSBarry Smith PetscFunctionBegin; 11594482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 116070441072SBarry Smith if (abstol) *abstol = snes->abstol; 116133174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 116233174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 116333174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 116433174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 11653a40ed3dSBarry Smith PetscFunctionReturn(0); 116633174efeSLois Curfman McInnes } 116733174efeSLois Curfman McInnes 11684a2ae208SSatish Balay #undef __FUNCT__ 11694a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 117033174efeSLois Curfman McInnes /*@ 11719b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 11729b94acceSBarry Smith 1173fee21e36SBarry Smith Collective on SNES 1174fee21e36SBarry Smith 1175c7afd0dbSLois Curfman McInnes Input Parameters: 1176c7afd0dbSLois Curfman McInnes + snes - the SNES context 1177c7afd0dbSLois Curfman McInnes - tol - tolerance 1178c7afd0dbSLois Curfman McInnes 11799b94acceSBarry Smith Options Database Key: 1180c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 11819b94acceSBarry Smith 118236851e7fSLois Curfman McInnes Level: intermediate 118336851e7fSLois Curfman McInnes 11849b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 11859b94acceSBarry Smith 11862492ecdbSBarry Smith .seealso: SNESSetTolerances() 11879b94acceSBarry Smith @*/ 118863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 11899b94acceSBarry Smith { 11903a40ed3dSBarry Smith PetscFunctionBegin; 11914482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 11929b94acceSBarry Smith snes->deltatol = tol; 11933a40ed3dSBarry Smith PetscFunctionReturn(0); 11949b94acceSBarry Smith } 11959b94acceSBarry Smith 1196df9fa365SBarry Smith /* 1197df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1198df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1199df9fa365SBarry Smith macros instead of functions 1200df9fa365SBarry Smith */ 12014a2ae208SSatish Balay #undef __FUNCT__ 12024a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 120363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1204ce1608b8SBarry Smith { 1205dfbe8321SBarry Smith PetscErrorCode ierr; 1206ce1608b8SBarry Smith 1207ce1608b8SBarry Smith PetscFunctionBegin; 12084482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1209ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1210ce1608b8SBarry Smith PetscFunctionReturn(0); 1211ce1608b8SBarry Smith } 1212ce1608b8SBarry Smith 12134a2ae208SSatish Balay #undef __FUNCT__ 12144a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 121563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1216df9fa365SBarry Smith { 1217dfbe8321SBarry Smith PetscErrorCode ierr; 1218df9fa365SBarry Smith 1219df9fa365SBarry Smith PetscFunctionBegin; 1220df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1221df9fa365SBarry Smith PetscFunctionReturn(0); 1222df9fa365SBarry Smith } 1223df9fa365SBarry Smith 12244a2ae208SSatish Balay #undef __FUNCT__ 12254a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 122663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw) 1227df9fa365SBarry Smith { 1228dfbe8321SBarry Smith PetscErrorCode ierr; 1229df9fa365SBarry Smith 1230df9fa365SBarry Smith PetscFunctionBegin; 1231df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1232df9fa365SBarry Smith PetscFunctionReturn(0); 1233df9fa365SBarry Smith } 1234df9fa365SBarry Smith 12359b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12369b94acceSBarry Smith 12374a2ae208SSatish Balay #undef __FUNCT__ 12384a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12399b94acceSBarry Smith /*@C 12405cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12419b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12429b94acceSBarry Smith progress. 12439b94acceSBarry Smith 1244fee21e36SBarry Smith Collective on SNES 1245fee21e36SBarry Smith 1246c7afd0dbSLois Curfman McInnes Input Parameters: 1247c7afd0dbSLois Curfman McInnes + snes - the SNES context 1248c7afd0dbSLois Curfman McInnes . func - monitoring routine 1249b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1250b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desitre) 1251b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1252b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 12539b94acceSBarry Smith 1254c7afd0dbSLois Curfman McInnes Calling sequence of func: 1255a7cc72afSBarry Smith $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1256c7afd0dbSLois Curfman McInnes 1257c7afd0dbSLois Curfman McInnes + snes - the SNES context 1258c7afd0dbSLois Curfman McInnes . its - iteration number 1259c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 126040a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 12619b94acceSBarry Smith 12629665c990SLois Curfman McInnes Options Database Keys: 1263c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1264c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1265c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1266c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1267c7afd0dbSLois Curfman McInnes been hardwired into a code by 1268c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1269c7afd0dbSLois Curfman McInnes does not cancel those set via 1270c7afd0dbSLois Curfman McInnes the options database. 12719665c990SLois Curfman McInnes 1272639f9d9dSBarry Smith Notes: 12736bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 12746bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 12756bc08f3fSLois Curfman McInnes order in which they were set. 1276639f9d9dSBarry Smith 127736851e7fSLois Curfman McInnes Level: intermediate 127836851e7fSLois Curfman McInnes 12799b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 12809b94acceSBarry Smith 12815cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 12829b94acceSBarry Smith @*/ 128363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 12849b94acceSBarry Smith { 12853a40ed3dSBarry Smith PetscFunctionBegin; 12864482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1287639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 128829bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1289639f9d9dSBarry Smith } 1290639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1291b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1292639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 12933a40ed3dSBarry Smith PetscFunctionReturn(0); 12949b94acceSBarry Smith } 12959b94acceSBarry Smith 12964a2ae208SSatish Balay #undef __FUNCT__ 12974a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 12985cd90555SBarry Smith /*@C 12995cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 13005cd90555SBarry Smith 1301c7afd0dbSLois Curfman McInnes Collective on SNES 1302c7afd0dbSLois Curfman McInnes 13035cd90555SBarry Smith Input Parameters: 13045cd90555SBarry Smith . snes - the SNES context 13055cd90555SBarry Smith 13061a480d89SAdministrator Options Database Key: 1307c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1308c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1309c7afd0dbSLois Curfman McInnes set via the options database 13105cd90555SBarry Smith 13115cd90555SBarry Smith Notes: 13125cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 13135cd90555SBarry Smith 131436851e7fSLois Curfman McInnes Level: intermediate 131536851e7fSLois Curfman McInnes 13165cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 13175cd90555SBarry Smith 13185cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 13195cd90555SBarry Smith @*/ 132063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes) 13215cd90555SBarry Smith { 13225cd90555SBarry Smith PetscFunctionBegin; 13234482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13245cd90555SBarry Smith snes->numbermonitors = 0; 13255cd90555SBarry Smith PetscFunctionReturn(0); 13265cd90555SBarry Smith } 13275cd90555SBarry Smith 13284a2ae208SSatish Balay #undef __FUNCT__ 13294a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13309b94acceSBarry Smith /*@C 13319b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13329b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13339b94acceSBarry Smith 1334fee21e36SBarry Smith Collective on SNES 1335fee21e36SBarry Smith 1336c7afd0dbSLois Curfman McInnes Input Parameters: 1337c7afd0dbSLois Curfman McInnes + snes - the SNES context 1338c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1339c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1340c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 13419b94acceSBarry Smith 1342c7afd0dbSLois Curfman McInnes Calling sequence of func: 13436849ba73SBarry Smith $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1344c7afd0dbSLois Curfman McInnes 1345c7afd0dbSLois Curfman McInnes + snes - the SNES context 1346c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1347184914b5SBarry Smith . reason - reason for convergence/divergence 1348c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 13494b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 13504b27c08aSLois Curfman McInnes - f - 2-norm of function 13519b94acceSBarry Smith 135236851e7fSLois Curfman McInnes Level: advanced 135336851e7fSLois Curfman McInnes 13549b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13559b94acceSBarry Smith 13564b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 13579b94acceSBarry Smith @*/ 135863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 13599b94acceSBarry Smith { 13603a40ed3dSBarry Smith PetscFunctionBegin; 13614482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13629b94acceSBarry Smith (snes)->converged = func; 13639b94acceSBarry Smith (snes)->cnvP = cctx; 13643a40ed3dSBarry Smith PetscFunctionReturn(0); 13659b94acceSBarry Smith } 13669b94acceSBarry Smith 13674a2ae208SSatish Balay #undef __FUNCT__ 13684a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 1369184914b5SBarry Smith /*@C 1370184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1371184914b5SBarry Smith 1372184914b5SBarry Smith Not Collective 1373184914b5SBarry Smith 1374184914b5SBarry Smith Input Parameter: 1375184914b5SBarry Smith . snes - the SNES context 1376184914b5SBarry Smith 1377184914b5SBarry Smith Output Parameter: 1378e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1379184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1380184914b5SBarry Smith 1381184914b5SBarry Smith Level: intermediate 1382184914b5SBarry Smith 1383184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1384184914b5SBarry Smith 1385184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1386184914b5SBarry Smith 13874b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1388184914b5SBarry Smith @*/ 138963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1390184914b5SBarry Smith { 1391184914b5SBarry Smith PetscFunctionBegin; 13924482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13934482741eSBarry Smith PetscValidPointer(reason,2); 1394184914b5SBarry Smith *reason = snes->reason; 1395184914b5SBarry Smith PetscFunctionReturn(0); 1396184914b5SBarry Smith } 1397184914b5SBarry Smith 13984a2ae208SSatish Balay #undef __FUNCT__ 13994a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1400c9005455SLois Curfman McInnes /*@ 1401c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1402c9005455SLois Curfman McInnes 1403fee21e36SBarry Smith Collective on SNES 1404fee21e36SBarry Smith 1405c7afd0dbSLois Curfman McInnes Input Parameters: 1406c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1407c7afd0dbSLois Curfman McInnes . a - array to hold history 1408cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1409758f92a0SBarry Smith . na - size of a and its 141064731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1411758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1412c7afd0dbSLois Curfman McInnes 1413c9005455SLois Curfman McInnes Notes: 14144b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1415c9005455SLois Curfman McInnes at each step. 1416c9005455SLois Curfman McInnes 1417c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1418c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1419c9005455SLois Curfman McInnes during the section of code that is being timed. 1420c9005455SLois Curfman McInnes 142136851e7fSLois Curfman McInnes Level: intermediate 142236851e7fSLois Curfman McInnes 1423c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1424758f92a0SBarry Smith 142508405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1426758f92a0SBarry Smith 1427c9005455SLois Curfman McInnes @*/ 142863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1429c9005455SLois Curfman McInnes { 14303a40ed3dSBarry Smith PetscFunctionBegin; 14314482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14324482741eSBarry Smith if (na) PetscValidScalarPointer(a,2); 1433c9005455SLois Curfman McInnes snes->conv_hist = a; 1434758f92a0SBarry Smith snes->conv_hist_its = its; 1435758f92a0SBarry Smith snes->conv_hist_max = na; 1436758f92a0SBarry Smith snes->conv_hist_reset = reset; 1437758f92a0SBarry Smith PetscFunctionReturn(0); 1438758f92a0SBarry Smith } 1439758f92a0SBarry Smith 14404a2ae208SSatish Balay #undef __FUNCT__ 14414a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 14420c4c9dddSBarry Smith /*@C 1443758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1444758f92a0SBarry Smith 1445758f92a0SBarry Smith Collective on SNES 1446758f92a0SBarry Smith 1447758f92a0SBarry Smith Input Parameter: 1448758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1449758f92a0SBarry Smith 1450758f92a0SBarry Smith Output Parameters: 1451758f92a0SBarry Smith . a - array to hold history 1452758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1453758f92a0SBarry Smith negative if not converged) for each solve. 1454758f92a0SBarry Smith - na - size of a and its 1455758f92a0SBarry Smith 1456758f92a0SBarry Smith Notes: 1457758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1458758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1459758f92a0SBarry Smith 1460758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1461758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1462758f92a0SBarry Smith during the section of code that is being timed. 1463758f92a0SBarry Smith 1464758f92a0SBarry Smith Level: intermediate 1465758f92a0SBarry Smith 1466758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1467758f92a0SBarry Smith 1468758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1469758f92a0SBarry Smith 1470758f92a0SBarry Smith @*/ 147163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1472758f92a0SBarry Smith { 1473758f92a0SBarry Smith PetscFunctionBegin; 14744482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1475758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1476758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1477758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 14783a40ed3dSBarry Smith PetscFunctionReturn(0); 1479c9005455SLois Curfman McInnes } 1480c9005455SLois Curfman McInnes 1481e74ef692SMatthew Knepley #undef __FUNCT__ 1482e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 1483ac226902SBarry Smith /*@C 148476b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 148576b2cf59SMatthew Knepley at the beginning of every step of the iteration. 148676b2cf59SMatthew Knepley 148776b2cf59SMatthew Knepley Collective on SNES 148876b2cf59SMatthew Knepley 148976b2cf59SMatthew Knepley Input Parameters: 149076b2cf59SMatthew Knepley . snes - The nonlinear solver context 149176b2cf59SMatthew Knepley . func - The function 149276b2cf59SMatthew Knepley 149376b2cf59SMatthew Knepley Calling sequence of func: 1494b5d30489SBarry Smith . func (SNES snes, PetscInt step); 149576b2cf59SMatthew Knepley 149676b2cf59SMatthew Knepley . step - The current step of the iteration 149776b2cf59SMatthew Knepley 149876b2cf59SMatthew Knepley Level: intermediate 149976b2cf59SMatthew Knepley 150076b2cf59SMatthew Knepley .keywords: SNES, update 1501b5d30489SBarry Smith 150276b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 150376b2cf59SMatthew Knepley @*/ 150463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 150576b2cf59SMatthew Knepley { 150676b2cf59SMatthew Knepley PetscFunctionBegin; 15074482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 150876b2cf59SMatthew Knepley snes->update = func; 150976b2cf59SMatthew Knepley PetscFunctionReturn(0); 151076b2cf59SMatthew Knepley } 151176b2cf59SMatthew Knepley 1512e74ef692SMatthew Knepley #undef __FUNCT__ 1513e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 151476b2cf59SMatthew Knepley /*@ 151576b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 151676b2cf59SMatthew Knepley 151776b2cf59SMatthew Knepley Not collective 151876b2cf59SMatthew Knepley 151976b2cf59SMatthew Knepley Input Parameters: 152076b2cf59SMatthew Knepley . snes - The nonlinear solver context 152176b2cf59SMatthew Knepley . step - The current step of the iteration 152276b2cf59SMatthew Knepley 1523205452f4SMatthew Knepley Level: intermediate 1524205452f4SMatthew Knepley 152576b2cf59SMatthew Knepley .keywords: SNES, update 152676b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 152776b2cf59SMatthew Knepley @*/ 152863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 152976b2cf59SMatthew Knepley { 153076b2cf59SMatthew Knepley PetscFunctionBegin; 153176b2cf59SMatthew Knepley PetscFunctionReturn(0); 153276b2cf59SMatthew Knepley } 153376b2cf59SMatthew Knepley 15344a2ae208SSatish Balay #undef __FUNCT__ 15354a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 15369b94acceSBarry Smith /* 15379b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 15389b94acceSBarry Smith positive parameter delta. 15399b94acceSBarry Smith 15409b94acceSBarry Smith Input Parameters: 1541c7afd0dbSLois Curfman McInnes + snes - the SNES context 15429b94acceSBarry Smith . y - approximate solution of linear system 15439b94acceSBarry Smith . fnorm - 2-norm of current function 1544c7afd0dbSLois Curfman McInnes - delta - trust region size 15459b94acceSBarry Smith 15469b94acceSBarry Smith Output Parameters: 1547c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 15489b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 15499b94acceSBarry Smith region, and exceeds zero otherwise. 1550c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 15519b94acceSBarry Smith 15529b94acceSBarry Smith Note: 15534b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 15549b94acceSBarry Smith is set to be the maximum allowable step size. 15559b94acceSBarry Smith 15569b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 15579b94acceSBarry Smith */ 1558dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 15599b94acceSBarry Smith { 1560064f8208SBarry Smith PetscReal nrm; 1561ea709b57SSatish Balay PetscScalar cnorm; 1562dfbe8321SBarry Smith PetscErrorCode ierr; 15633a40ed3dSBarry Smith 15643a40ed3dSBarry Smith PetscFunctionBegin; 15654482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 15664482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1567c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,2); 1568184914b5SBarry Smith 1569064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1570064f8208SBarry Smith if (nrm > *delta) { 1571064f8208SBarry Smith nrm = *delta/nrm; 1572064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1573064f8208SBarry Smith cnorm = nrm; 1574329f5518SBarry Smith ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 15759b94acceSBarry Smith *ynorm = *delta; 15769b94acceSBarry Smith } else { 15779b94acceSBarry Smith *gpnorm = 0.0; 1578064f8208SBarry Smith *ynorm = nrm; 15799b94acceSBarry Smith } 15803a40ed3dSBarry Smith PetscFunctionReturn(0); 15819b94acceSBarry Smith } 15829b94acceSBarry Smith 15834a2ae208SSatish Balay #undef __FUNCT__ 15844a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 15859b94acceSBarry Smith /*@ 15869b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 158763a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 15889b94acceSBarry Smith 1589c7afd0dbSLois Curfman McInnes Collective on SNES 1590c7afd0dbSLois Curfman McInnes 1591b2002411SLois Curfman McInnes Input Parameters: 1592c7afd0dbSLois Curfman McInnes + snes - the SNES context 159370e92668SMatthew Knepley - x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution() 15949b94acceSBarry Smith 1595b2002411SLois Curfman McInnes Notes: 15968ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 15978ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 15988ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 15998ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 16008ddd3da0SLois Curfman McInnes 160136851e7fSLois Curfman McInnes Level: beginner 160236851e7fSLois Curfman McInnes 16039b94acceSBarry Smith .keywords: SNES, nonlinear, solve 16049b94acceSBarry Smith 160570e92668SMatthew Knepley .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution() 16069b94acceSBarry Smith @*/ 160763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec x) 16089b94acceSBarry Smith { 1609dfbe8321SBarry Smith PetscErrorCode ierr; 1610f1af5d2fSBarry Smith PetscTruth flg; 1611052efed2SBarry Smith 16123a40ed3dSBarry Smith PetscFunctionBegin; 16134482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16141302d50aSBarry Smith if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 161574637425SBarry Smith 161670e92668SMatthew Knepley if (x) { 161770e92668SMatthew Knepley PetscValidHeaderSpecific(x,VEC_COOKIE,2); 161870e92668SMatthew Knepley PetscCheckSameComm(snes,1,x,2); 161970e92668SMatthew Knepley } else { 162070e92668SMatthew Knepley ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr); 162170e92668SMatthew Knepley if (!x) { 162270e92668SMatthew Knepley ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr); 162370e92668SMatthew Knepley } 162470e92668SMatthew Knepley } 162570e92668SMatthew Knepley snes->vec_sol = snes->vec_sol_always = x; 162670e92668SMatthew Knepley if (!snes->setupcalled) { 162770e92668SMatthew Knepley ierr = SNESSetUp(snes);CHKERRQ(ierr); 162870e92668SMatthew Knepley } 1629abc0a331SBarry Smith if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1630d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 163150ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1632d5e45103SBarry Smith 1633d5e45103SBarry Smith ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN); 1634d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 163538f152feSBarry Smith /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */ 163619717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr); 1637d5e45103SBarry Smith } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) { 1638d5e45103SBarry Smith /* translate exception into SNES not converged reason */ 1639d5e45103SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 164038f152feSBarry Smith ierr = 0; 164138f152feSBarry Smith } 164238f152feSBarry Smith CHKERRQ(ierr); 1643d5e45103SBarry Smith 1644d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1645b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 16468b179ff8SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1647da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1648da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 16495968eb51SBarry Smith if (snes->printreason) { 16505968eb51SBarry Smith if (snes->reason > 0) { 16519dcbbd2bSBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 16525968eb51SBarry Smith } else { 16539dcbbd2bSBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 16545968eb51SBarry Smith } 16555968eb51SBarry Smith } 16565968eb51SBarry Smith 16573a40ed3dSBarry Smith PetscFunctionReturn(0); 16589b94acceSBarry Smith } 16599b94acceSBarry Smith 16609b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 16619b94acceSBarry Smith 16624a2ae208SSatish Balay #undef __FUNCT__ 16634a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 166482bf6240SBarry Smith /*@C 16654b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 16669b94acceSBarry Smith 1667fee21e36SBarry Smith Collective on SNES 1668fee21e36SBarry Smith 1669c7afd0dbSLois Curfman McInnes Input Parameters: 1670c7afd0dbSLois Curfman McInnes + snes - the SNES context 1671454a90a3SBarry Smith - type - a known method 1672c7afd0dbSLois Curfman McInnes 1673c7afd0dbSLois Curfman McInnes Options Database Key: 1674454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1675c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1676ae12b187SLois Curfman McInnes 16779b94acceSBarry Smith Notes: 1678e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 16794b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1680c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 16814b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1682c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 16839b94acceSBarry Smith 1684ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1685ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1686ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1687ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1688ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1689ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1690ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1691ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1692ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1693b0a32e0cSBarry Smith appropriate method. 169436851e7fSLois Curfman McInnes 169536851e7fSLois Curfman McInnes Level: intermediate 1696a703fe33SLois Curfman McInnes 1697454a90a3SBarry Smith .keywords: SNES, set, type 1698435da068SBarry Smith 1699435da068SBarry Smith .seealso: SNESType, SNESCreate() 1700435da068SBarry Smith 17019b94acceSBarry Smith @*/ 170263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type) 17039b94acceSBarry Smith { 1704dfbe8321SBarry Smith PetscErrorCode ierr,(*r)(SNES); 17056831982aSBarry Smith PetscTruth match; 17063a40ed3dSBarry Smith 17073a40ed3dSBarry Smith PetscFunctionBegin; 17084482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17094482741eSBarry Smith PetscValidCharPointer(type,2); 171082bf6240SBarry Smith 17116831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 17120f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 171382bf6240SBarry Smith 171482bf6240SBarry Smith if (snes->setupcalled) { 1715e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 171682bf6240SBarry Smith snes->data = 0; 171782bf6240SBarry Smith } 17187d1a2b2bSBarry Smith 17199b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 172082bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1721b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1722958c9bccSBarry Smith if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1723606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 172482bf6240SBarry Smith snes->data = 0; 17253a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1726454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 17273a40ed3dSBarry Smith PetscFunctionReturn(0); 17289b94acceSBarry Smith } 17299b94acceSBarry Smith 1730a847f771SSatish Balay 17319b94acceSBarry Smith /* --------------------------------------------------------------------- */ 17324a2ae208SSatish Balay #undef __FUNCT__ 17334a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 17349b94acceSBarry Smith /*@C 17359b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1736f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 17379b94acceSBarry Smith 1738fee21e36SBarry Smith Not Collective 1739fee21e36SBarry Smith 174036851e7fSLois Curfman McInnes Level: advanced 174136851e7fSLois Curfman McInnes 17429b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 17439b94acceSBarry Smith 17449b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 17459b94acceSBarry Smith @*/ 174663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 17479b94acceSBarry Smith { 1748dfbe8321SBarry Smith PetscErrorCode ierr; 174982bf6240SBarry Smith 17503a40ed3dSBarry Smith PetscFunctionBegin; 175182bf6240SBarry Smith if (SNESList) { 1752b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 175382bf6240SBarry Smith SNESList = 0; 17549b94acceSBarry Smith } 17554c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 17563a40ed3dSBarry Smith PetscFunctionReturn(0); 17579b94acceSBarry Smith } 17589b94acceSBarry Smith 17594a2ae208SSatish Balay #undef __FUNCT__ 17604a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 17619b94acceSBarry Smith /*@C 17629a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 17639b94acceSBarry Smith 1764c7afd0dbSLois Curfman McInnes Not Collective 1765c7afd0dbSLois Curfman McInnes 17669b94acceSBarry Smith Input Parameter: 17674b0e389bSBarry Smith . snes - nonlinear solver context 17689b94acceSBarry Smith 17699b94acceSBarry Smith Output Parameter: 17703a7fca6bSBarry Smith . type - SNES method (a character string) 17719b94acceSBarry Smith 177236851e7fSLois Curfman McInnes Level: intermediate 177336851e7fSLois Curfman McInnes 1774454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 17759b94acceSBarry Smith @*/ 177663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type) 17779b94acceSBarry Smith { 17783a40ed3dSBarry Smith PetscFunctionBegin; 17794482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17804482741eSBarry Smith PetscValidPointer(type,2); 1781454a90a3SBarry Smith *type = snes->type_name; 17823a40ed3dSBarry Smith PetscFunctionReturn(0); 17839b94acceSBarry Smith } 17849b94acceSBarry Smith 17854a2ae208SSatish Balay #undef __FUNCT__ 17864a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 17879b94acceSBarry Smith /*@C 17889b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 17899b94acceSBarry Smith stored. 17909b94acceSBarry Smith 1791c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1792c7afd0dbSLois Curfman McInnes 17939b94acceSBarry Smith Input Parameter: 17949b94acceSBarry Smith . snes - the SNES context 17959b94acceSBarry Smith 17969b94acceSBarry Smith Output Parameter: 17979b94acceSBarry Smith . x - the solution 17989b94acceSBarry Smith 179970e92668SMatthew Knepley Level: intermediate 180036851e7fSLois Curfman McInnes 18019b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 18029b94acceSBarry Smith 180370e92668SMatthew Knepley .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 18049b94acceSBarry Smith @*/ 180563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 18069b94acceSBarry Smith { 18073a40ed3dSBarry Smith PetscFunctionBegin; 18084482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18094482741eSBarry Smith PetscValidPointer(x,2); 18109b94acceSBarry Smith *x = snes->vec_sol_always; 18113a40ed3dSBarry Smith PetscFunctionReturn(0); 18129b94acceSBarry Smith } 18139b94acceSBarry Smith 18144a2ae208SSatish Balay #undef __FUNCT__ 181570e92668SMatthew Knepley #define __FUNCT__ "SNESSetSolution" 181670e92668SMatthew Knepley /*@C 181770e92668SMatthew Knepley SNESSetSolution - Sets the vector where the approximate solution is stored. 181870e92668SMatthew Knepley 181970e92668SMatthew Knepley Not Collective, but Vec is parallel if SNES is parallel 182070e92668SMatthew Knepley 182170e92668SMatthew Knepley Input Parameters: 182270e92668SMatthew Knepley + snes - the SNES context 182370e92668SMatthew Knepley - x - the solution 182470e92668SMatthew Knepley 182570e92668SMatthew Knepley Output Parameter: 182670e92668SMatthew Knepley 182770e92668SMatthew Knepley Level: intermediate 182870e92668SMatthew Knepley 182970e92668SMatthew Knepley .keywords: SNES, nonlinear, set, solution 183070e92668SMatthew Knepley 183170e92668SMatthew Knepley .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 183270e92668SMatthew Knepley @*/ 183370e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x) 183470e92668SMatthew Knepley { 183570e92668SMatthew Knepley PetscFunctionBegin; 183670e92668SMatthew Knepley PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 183770e92668SMatthew Knepley PetscValidHeaderSpecific(x,VEC_COOKIE,2); 183870e92668SMatthew Knepley PetscCheckSameComm(snes,1,x,2); 183970e92668SMatthew Knepley snes->vec_sol_always = x; 184070e92668SMatthew Knepley PetscFunctionReturn(0); 184170e92668SMatthew Knepley } 184270e92668SMatthew Knepley 184370e92668SMatthew Knepley #undef __FUNCT__ 18444a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 18459b94acceSBarry Smith /*@C 18469b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 18479b94acceSBarry Smith stored. 18489b94acceSBarry Smith 1849c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1850c7afd0dbSLois Curfman McInnes 18519b94acceSBarry Smith Input Parameter: 18529b94acceSBarry Smith . snes - the SNES context 18539b94acceSBarry Smith 18549b94acceSBarry Smith Output Parameter: 18559b94acceSBarry Smith . x - the solution update 18569b94acceSBarry Smith 185736851e7fSLois Curfman McInnes Level: advanced 185836851e7fSLois Curfman McInnes 18599b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 18609b94acceSBarry Smith 18619b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 18629b94acceSBarry Smith @*/ 186363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 18649b94acceSBarry Smith { 18653a40ed3dSBarry Smith PetscFunctionBegin; 18664482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18674482741eSBarry Smith PetscValidPointer(x,2); 18689b94acceSBarry Smith *x = snes->vec_sol_update_always; 18693a40ed3dSBarry Smith PetscFunctionReturn(0); 18709b94acceSBarry Smith } 18719b94acceSBarry Smith 18724a2ae208SSatish Balay #undef __FUNCT__ 18734a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 18749b94acceSBarry Smith /*@C 18753638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 18769b94acceSBarry Smith 1877c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1878c7afd0dbSLois Curfman McInnes 18799b94acceSBarry Smith Input Parameter: 18809b94acceSBarry Smith . snes - the SNES context 18819b94acceSBarry Smith 18829b94acceSBarry Smith Output Parameter: 18837bf4e008SBarry Smith + r - the function (or PETSC_NULL) 188470e92668SMatthew Knepley . func - the function (or PETSC_NULL) 188570e92668SMatthew Knepley - ctx - the function context (or PETSC_NULL) 18869b94acceSBarry Smith 188736851e7fSLois Curfman McInnes Level: advanced 188836851e7fSLois Curfman McInnes 1889a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 18909b94acceSBarry Smith 18914b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 18929b94acceSBarry Smith @*/ 189370e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 18949b94acceSBarry Smith { 18953a40ed3dSBarry Smith PetscFunctionBegin; 18964482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18977bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 189800036973SBarry Smith if (func) *func = snes->computefunction; 189970e92668SMatthew Knepley if (ctx) *ctx = snes->funP; 19003a40ed3dSBarry Smith PetscFunctionReturn(0); 19019b94acceSBarry Smith } 19029b94acceSBarry Smith 19034a2ae208SSatish Balay #undef __FUNCT__ 19044a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 19053c7409f5SSatish Balay /*@C 19063c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1907d850072dSLois Curfman McInnes SNES options in the database. 19083c7409f5SSatish Balay 1909fee21e36SBarry Smith Collective on SNES 1910fee21e36SBarry Smith 1911c7afd0dbSLois Curfman McInnes Input Parameter: 1912c7afd0dbSLois Curfman McInnes + snes - the SNES context 1913c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1914c7afd0dbSLois Curfman McInnes 1915d850072dSLois Curfman McInnes Notes: 1916a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1917c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1918d850072dSLois Curfman McInnes 191936851e7fSLois Curfman McInnes Level: advanced 192036851e7fSLois Curfman McInnes 19213c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1922a86d99e1SLois Curfman McInnes 1923a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19243c7409f5SSatish Balay @*/ 192563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 19263c7409f5SSatish Balay { 1927dfbe8321SBarry Smith PetscErrorCode ierr; 19283c7409f5SSatish Balay 19293a40ed3dSBarry Smith PetscFunctionBegin; 19304482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1931639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 193294b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 19333a40ed3dSBarry Smith PetscFunctionReturn(0); 19343c7409f5SSatish Balay } 19353c7409f5SSatish Balay 19364a2ae208SSatish Balay #undef __FUNCT__ 19374a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 19383c7409f5SSatish Balay /*@C 1939f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1940d850072dSLois Curfman McInnes SNES options in the database. 19413c7409f5SSatish Balay 1942fee21e36SBarry Smith Collective on SNES 1943fee21e36SBarry Smith 1944c7afd0dbSLois Curfman McInnes Input Parameters: 1945c7afd0dbSLois Curfman McInnes + snes - the SNES context 1946c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1947c7afd0dbSLois Curfman McInnes 1948d850072dSLois Curfman McInnes Notes: 1949a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1950c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1951d850072dSLois Curfman McInnes 195236851e7fSLois Curfman McInnes Level: advanced 195336851e7fSLois Curfman McInnes 19543c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1955a86d99e1SLois Curfman McInnes 1956a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 19573c7409f5SSatish Balay @*/ 195863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 19593c7409f5SSatish Balay { 1960dfbe8321SBarry Smith PetscErrorCode ierr; 19613c7409f5SSatish Balay 19623a40ed3dSBarry Smith PetscFunctionBegin; 19634482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1964639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 196594b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 19663a40ed3dSBarry Smith PetscFunctionReturn(0); 19673c7409f5SSatish Balay } 19683c7409f5SSatish Balay 19694a2ae208SSatish Balay #undef __FUNCT__ 19704a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 19719ab63eb5SSatish Balay /*@C 19723c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 19733c7409f5SSatish Balay SNES options in the database. 19743c7409f5SSatish Balay 1975c7afd0dbSLois Curfman McInnes Not Collective 1976c7afd0dbSLois Curfman McInnes 19773c7409f5SSatish Balay Input Parameter: 19783c7409f5SSatish Balay . snes - the SNES context 19793c7409f5SSatish Balay 19803c7409f5SSatish Balay Output Parameter: 19813c7409f5SSatish Balay . prefix - pointer to the prefix string used 19823c7409f5SSatish Balay 19839ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 19849ab63eb5SSatish Balay sufficient length to hold the prefix. 19859ab63eb5SSatish Balay 198636851e7fSLois Curfman McInnes Level: advanced 198736851e7fSLois Curfman McInnes 19883c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 1989a86d99e1SLois Curfman McInnes 1990a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 19913c7409f5SSatish Balay @*/ 199263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,char *prefix[]) 19933c7409f5SSatish Balay { 1994dfbe8321SBarry Smith PetscErrorCode ierr; 19953c7409f5SSatish Balay 19963a40ed3dSBarry Smith PetscFunctionBegin; 19974482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1998639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 19993a40ed3dSBarry Smith PetscFunctionReturn(0); 20003c7409f5SSatish Balay } 20013c7409f5SSatish Balay 2002b2002411SLois Curfman McInnes 20034a2ae208SSatish Balay #undef __FUNCT__ 20044a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 20053cea93caSBarry Smith /*@C 20063cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 20073cea93caSBarry Smith 20087f6c08e0SMatthew Knepley Level: advanced 20093cea93caSBarry Smith @*/ 201063dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2011b2002411SLois Curfman McInnes { 2012e2d1d2b7SBarry Smith char fullname[PETSC_MAX_PATH_LEN]; 2013dfbe8321SBarry Smith PetscErrorCode ierr; 2014b2002411SLois Curfman McInnes 2015b2002411SLois Curfman McInnes PetscFunctionBegin; 2016b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2017c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2018b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2019b2002411SLois Curfman McInnes } 2020da9b6338SBarry Smith 2021da9b6338SBarry Smith #undef __FUNCT__ 2022da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 202363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2024da9b6338SBarry Smith { 2025dfbe8321SBarry Smith PetscErrorCode ierr; 202677431f27SBarry Smith PetscInt N,i,j; 2027da9b6338SBarry Smith Vec u,uh,fh; 2028da9b6338SBarry Smith PetscScalar value; 2029da9b6338SBarry Smith PetscReal norm; 2030da9b6338SBarry Smith 2031da9b6338SBarry Smith PetscFunctionBegin; 2032da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2033da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2034da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2035da9b6338SBarry Smith 2036da9b6338SBarry Smith /* currently only works for sequential */ 2037da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2038da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2039da9b6338SBarry Smith for (i=0; i<N; i++) { 2040da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 204177431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2042da9b6338SBarry Smith for (j=-10; j<11; j++) { 2043ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2044da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 20453ab0aad5SBarry Smith ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2046da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 204777431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2048da9b6338SBarry Smith value = -value; 2049da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2050da9b6338SBarry Smith } 2051da9b6338SBarry Smith } 2052da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2053da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2054da9b6338SBarry Smith PetscFunctionReturn(0); 2055da9b6338SBarry Smith } 2056