163dd3a1aSKris Buschelman #define PETSCSNES_DLL 29b94acceSBarry Smith 3e090d566SSatish Balay #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 49b94acceSBarry Smith 54c49b128SBarry Smith PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 68ba1e511SMatthew Knepley PetscFList SNESList = PETSC_NULL; 78ba1e511SMatthew Knepley 88ba1e511SMatthew Knepley /* Logging support */ 963dd3a1aSKris Buschelman PetscCookie PETSCSNES_DLLEXPORT SNES_COOKIE = 0; 106849ba73SBarry Smith PetscEvent SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0; 11a09944afSBarry Smith 12a09944afSBarry Smith #undef __FUNCT__ 134a2ae208SSatish Balay #define __FUNCT__ "SNESView" 147e2c5f70SBarry Smith /*@C 159b94acceSBarry Smith SNESView - Prints the SNES data structure. 169b94acceSBarry Smith 174c49b128SBarry Smith Collective on SNES 18fee21e36SBarry Smith 19c7afd0dbSLois Curfman McInnes Input Parameters: 20c7afd0dbSLois Curfman McInnes + SNES - the SNES context 21c7afd0dbSLois Curfman McInnes - viewer - visualization context 22c7afd0dbSLois Curfman McInnes 239b94acceSBarry Smith Options Database Key: 24c8a8ba5cSLois Curfman McInnes . -snes_view - Calls SNESView() at end of SNESSolve() 259b94acceSBarry Smith 269b94acceSBarry Smith Notes: 279b94acceSBarry Smith The available visualization contexts include 28b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 29b0a32e0cSBarry Smith - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 30c8a8ba5cSLois Curfman McInnes output where only the first processor opens 31c8a8ba5cSLois Curfman McInnes the file. All other processors send their 32c8a8ba5cSLois Curfman McInnes data to the first processor to print. 339b94acceSBarry Smith 343e081fefSLois Curfman McInnes The user can open an alternative visualization context with 35b0a32e0cSBarry Smith PetscViewerASCIIOpen() - output to a specified file. 369b94acceSBarry Smith 3736851e7fSLois Curfman McInnes Level: beginner 3836851e7fSLois Curfman McInnes 399b94acceSBarry Smith .keywords: SNES, view 409b94acceSBarry Smith 41b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen() 429b94acceSBarry Smith @*/ 4363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer) 449b94acceSBarry Smith { 459b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 46dfbe8321SBarry Smith PetscErrorCode ierr; 4794b7f48cSBarry Smith KSP ksp; 48e060cb09SBarry Smith SNESType type; 4932077d6dSBarry Smith PetscTruth iascii,isstring; 509b94acceSBarry Smith 513a40ed3dSBarry Smith PetscFunctionBegin; 524482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 53b0a32e0cSBarry Smith if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm); 544482741eSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 55c9780b6fSBarry Smith PetscCheckSameComm(snes,1,viewer,2); 5674679c65SBarry Smith 5732077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 58b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 5932077d6dSBarry Smith if (iascii) { 603a7fca6bSBarry Smith if (snes->prefix) { 613a7fca6bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr); 623a7fca6bSBarry Smith } else { 63b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 643a7fca6bSBarry Smith } 65454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 66454a90a3SBarry Smith if (type) { 67b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 68184914b5SBarry Smith } else { 69b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 70184914b5SBarry Smith } 710ef38995SBarry Smith if (snes->view) { 72b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 730ef38995SBarry Smith ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 74b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 750ef38995SBarry Smith } 7677431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 77a83599f4SBarry 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); 85a83599f4SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 86a83599f4SBarry 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 156e8105e01SRichard Katz . -snes_monitor <optional filename> - prints residual norm at each iteration. if no 157e8105e01SRichard Katz filename given prints to stdout 158d132466eSBarry Smith . -snes_vecmonitor - plots solution at each iteration 159d132466eSBarry Smith . -snes_vecmonitor_update - plots update to solution at each iteration 16082738288SBarry Smith . -snes_xmonitor - plots residual norm at each iteration 161e24b481bSBarry Smith . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 1625968eb51SBarry Smith . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 1635968eb51SBarry Smith - -snes_print_converged_reason - print the reason for convergence/divergence after each solve 16482738288SBarry Smith 16582738288SBarry Smith Options Database for Eisenstat-Walker method: 1664b27c08aSLois Curfman McInnes + -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence 1674b27c08aSLois Curfman McInnes . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 16836851e7fSLois Curfman McInnes . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 16936851e7fSLois Curfman McInnes . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 17036851e7fSLois Curfman McInnes . -snes_ksp_ew_gamma <gamma> - Sets gamma 17136851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha <alpha> - Sets alpha 17236851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 17336851e7fSLois Curfman McInnes - -snes_ksp_ew_threshold <threshold> - Sets threshold 17482738288SBarry Smith 17511ca99fdSLois Curfman McInnes Notes: 17611ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 17711ca99fdSLois Curfman McInnes the users manual. 17883e2fdc7SBarry Smith 17936851e7fSLois Curfman McInnes Level: beginner 18036851e7fSLois Curfman McInnes 1819b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1829b94acceSBarry Smith 18369ed319cSSatish Balay .seealso: SNESSetOptionsPrefix() 1849b94acceSBarry Smith @*/ 18563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes) 1869b94acceSBarry Smith { 18794b7f48cSBarry Smith KSP ksp; 188186905e3SBarry Smith SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 189f1af5d2fSBarry Smith PetscTruth flg; 190dfbe8321SBarry Smith PetscErrorCode ierr; 19177431f27SBarry Smith PetscInt i; 1922fc52814SBarry Smith const char *deft; 193e8105e01SRichard Katz char type[256], monfilename[PETSC_MAX_PATH_LEN]; 194e8105e01SRichard Katz PetscViewer monviewer; 1959b94acceSBarry Smith 1963a40ed3dSBarry Smith PetscFunctionBegin; 1974482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 198ca161407SBarry Smith 199b0a32e0cSBarry Smith ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 200186905e3SBarry Smith if (snes->type_name) { 201186905e3SBarry Smith deft = snes->type_name; 202186905e3SBarry Smith } else { 2034b27c08aSLois Curfman McInnes deft = SNESLS; 204d64ed03dSBarry Smith } 2054bbc92c1SBarry Smith 206186905e3SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 207b0a32e0cSBarry Smith ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 208d64ed03dSBarry Smith if (flg) { 209186905e3SBarry Smith ierr = SNESSetType(snes,type);CHKERRQ(ierr); 210186905e3SBarry Smith } else if (!snes->type_name) { 211186905e3SBarry Smith ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 212d64ed03dSBarry Smith } 213909c8a9fSBarry Smith ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 21493c39befSBarry Smith 21587828ca2SBarry Smith ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 21670441072SBarry Smith ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr); 217186905e3SBarry Smith 21887828ca2SBarry Smith ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 219b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 220b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 22150ffb88aSMatthew Knepley ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 2225968eb51SBarry Smith ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr); 2235968eb51SBarry Smith if (flg) { 2245968eb51SBarry Smith snes->printreason = PETSC_TRUE; 2255968eb51SBarry Smith } 226186905e3SBarry Smith 2277aaed0d8SBarry Smith ierr = PetscOptionsTruth("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr); 228186905e3SBarry Smith 229b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 23087828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 23187828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 23287828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 23387828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 23487828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 23587828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 236186905e3SBarry Smith 237b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr); 23893c39befSBarry Smith if (flg) {snes->converged = 0;} 239b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr); 2405cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 241eabae89aSBarry Smith 242eabae89aSBarry Smith ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 243e8105e01SRichard Katz if (flg) { 244e8105e01SRichard Katz ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr); 245e8105e01SRichard Katz ierr = SNESSetMonitor(snes,SNESDefaultMonitor,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 246e8105e01SRichard Katz } 247eabae89aSBarry Smith 248eabae89aSBarry Smith ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESSetRatioMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 249eabae89aSBarry Smith if (flg) { 250eabae89aSBarry Smith ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr); 251eabae89aSBarry Smith ierr = SNESSetRatioMonitor(snes,monviewer);CHKERRQ(ierr); 252e8105e01SRichard Katz } 253eabae89aSBarry Smith 254eabae89aSBarry Smith ierr = PetscOptionsString("-snes_smonitor","Monitor norm of function (fewer digits)","SNESSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 255eabae89aSBarry Smith if (flg) { 256eabae89aSBarry Smith ierr = PetscViewerASCIIOpen(snes->comm,monfilename,&monviewer);CHKERRQ(ierr); 257eabae89aSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,monviewer,(PetscErrorCode (*)(void*))PetscViewerDestroy);CHKERRQ(ierr); 258eabae89aSBarry Smith } 259eabae89aSBarry Smith 260b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr); 261b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 262b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr); 2637c922b88SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);} 2645ed2d596SBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr); 2655ed2d596SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);} 266b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr); 267186905e3SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 268e24b481bSBarry Smith 269b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 2704b27c08aSLois Curfman McInnes if (flg) { 271186905e3SBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 272ae15b995SBarry Smith ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr); 2739b94acceSBarry Smith } 274639f9d9dSBarry Smith 27576b2cf59SMatthew Knepley for(i = 0; i < numberofsetfromoptions; i++) { 27676b2cf59SMatthew Knepley ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 27776b2cf59SMatthew Knepley } 27876b2cf59SMatthew Knepley 279186905e3SBarry Smith if (snes->setfromoptions) { 280186905e3SBarry Smith ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 281639f9d9dSBarry Smith } 282b0a32e0cSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 2834bbc92c1SBarry Smith 28494b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 28594b7f48cSBarry Smith ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 28693993e2dSLois Curfman McInnes 2873a40ed3dSBarry Smith PetscFunctionReturn(0); 2889b94acceSBarry Smith } 2899b94acceSBarry Smith 290a847f771SSatish Balay 2914a2ae208SSatish Balay #undef __FUNCT__ 2924a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext" 2939b94acceSBarry Smith /*@ 2949b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 2959b94acceSBarry Smith the nonlinear solvers. 2969b94acceSBarry Smith 297fee21e36SBarry Smith Collective on SNES 298fee21e36SBarry Smith 299c7afd0dbSLois Curfman McInnes Input Parameters: 300c7afd0dbSLois Curfman McInnes + snes - the SNES context 301c7afd0dbSLois Curfman McInnes - usrP - optional user context 302c7afd0dbSLois Curfman McInnes 30336851e7fSLois Curfman McInnes Level: intermediate 30436851e7fSLois Curfman McInnes 3059b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3069b94acceSBarry Smith 3079b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3089b94acceSBarry Smith @*/ 30963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP) 3109b94acceSBarry Smith { 3113a40ed3dSBarry Smith PetscFunctionBegin; 3124482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 3139b94acceSBarry Smith snes->user = usrP; 3143a40ed3dSBarry Smith PetscFunctionReturn(0); 3159b94acceSBarry Smith } 31674679c65SBarry Smith 3174a2ae208SSatish Balay #undef __FUNCT__ 3184a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext" 3199b94acceSBarry Smith /*@C 3209b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3219b94acceSBarry Smith nonlinear solvers. 3229b94acceSBarry Smith 323c7afd0dbSLois Curfman McInnes Not Collective 324c7afd0dbSLois Curfman McInnes 3259b94acceSBarry Smith Input Parameter: 3269b94acceSBarry Smith . snes - SNES context 3279b94acceSBarry Smith 3289b94acceSBarry Smith Output Parameter: 3299b94acceSBarry Smith . usrP - user context 3309b94acceSBarry Smith 33136851e7fSLois Curfman McInnes Level: intermediate 33236851e7fSLois Curfman McInnes 3339b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3349b94acceSBarry Smith 3359b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3369b94acceSBarry Smith @*/ 33763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP) 3389b94acceSBarry Smith { 3393a40ed3dSBarry Smith PetscFunctionBegin; 3404482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 3419b94acceSBarry Smith *usrP = snes->user; 3423a40ed3dSBarry Smith PetscFunctionReturn(0); 3439b94acceSBarry Smith } 34474679c65SBarry Smith 3454a2ae208SSatish Balay #undef __FUNCT__ 3464a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber" 3479b94acceSBarry Smith /*@ 348c8228a4eSBarry Smith SNESGetIterationNumber - Gets the number of nonlinear iterations completed 349c8228a4eSBarry Smith at this time. 3509b94acceSBarry Smith 351c7afd0dbSLois Curfman McInnes Not Collective 352c7afd0dbSLois Curfman McInnes 3539b94acceSBarry Smith Input Parameter: 3549b94acceSBarry Smith . snes - SNES context 3559b94acceSBarry Smith 3569b94acceSBarry Smith Output Parameter: 3579b94acceSBarry Smith . iter - iteration number 3589b94acceSBarry Smith 359c8228a4eSBarry Smith Notes: 360c8228a4eSBarry Smith For example, during the computation of iteration 2 this would return 1. 361c8228a4eSBarry Smith 362c8228a4eSBarry Smith This is useful for using lagged Jacobians (where one does not recompute the 36308405cd6SLois Curfman McInnes Jacobian at each SNES iteration). For example, the code 36408405cd6SLois Curfman McInnes .vb 36508405cd6SLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&it); 36608405cd6SLois Curfman McInnes if (!(it % 2)) { 36708405cd6SLois Curfman McInnes [compute Jacobian here] 36808405cd6SLois Curfman McInnes } 36908405cd6SLois Curfman McInnes .ve 370c8228a4eSBarry Smith can be used in your ComputeJacobian() function to cause the Jacobian to be 37108405cd6SLois Curfman McInnes recomputed every second SNES iteration. 372c8228a4eSBarry Smith 37336851e7fSLois Curfman McInnes Level: intermediate 37436851e7fSLois Curfman McInnes 375*2b668275SBarry Smith .keywords: SNES, nonlinear, get, iteration, number, 376*2b668275SBarry Smith 377*2b668275SBarry Smith .seealso: SNESGetFunctionNorm(), SNESGetNumberLinearIterations() 3789b94acceSBarry Smith @*/ 37963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter) 3809b94acceSBarry Smith { 3813a40ed3dSBarry Smith PetscFunctionBegin; 3824482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 3834482741eSBarry Smith PetscValidIntPointer(iter,2); 3849b94acceSBarry Smith *iter = snes->iter; 3853a40ed3dSBarry Smith PetscFunctionReturn(0); 3869b94acceSBarry Smith } 38774679c65SBarry Smith 3884a2ae208SSatish Balay #undef __FUNCT__ 3894a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm" 3909b94acceSBarry Smith /*@ 3919b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3929b94acceSBarry Smith with SNESSSetFunction(). 3939b94acceSBarry Smith 394c7afd0dbSLois Curfman McInnes Collective on SNES 395c7afd0dbSLois Curfman McInnes 3969b94acceSBarry Smith Input Parameter: 3979b94acceSBarry Smith . snes - SNES context 3989b94acceSBarry Smith 3999b94acceSBarry Smith Output Parameter: 4009b94acceSBarry Smith . fnorm - 2-norm of function 4019b94acceSBarry Smith 40236851e7fSLois Curfman McInnes Level: intermediate 40336851e7fSLois Curfman McInnes 4049b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 405a86d99e1SLois Curfman McInnes 406*2b668275SBarry Smith .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetNumberLinearIterations() 4079b94acceSBarry Smith @*/ 40863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm) 4099b94acceSBarry Smith { 4103a40ed3dSBarry Smith PetscFunctionBegin; 4114482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 4124482741eSBarry Smith PetscValidScalarPointer(fnorm,2); 4139b94acceSBarry Smith *fnorm = snes->norm; 4143a40ed3dSBarry Smith PetscFunctionReturn(0); 4159b94acceSBarry Smith } 41674679c65SBarry Smith 4174a2ae208SSatish Balay #undef __FUNCT__ 4184a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps" 4199b94acceSBarry Smith /*@ 4209b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4219b94acceSBarry Smith attempted by the nonlinear solver. 4229b94acceSBarry Smith 423c7afd0dbSLois Curfman McInnes Not Collective 424c7afd0dbSLois Curfman McInnes 4259b94acceSBarry Smith Input Parameter: 4269b94acceSBarry Smith . snes - SNES context 4279b94acceSBarry Smith 4289b94acceSBarry Smith Output Parameter: 4299b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4309b94acceSBarry Smith 431c96a6f78SLois Curfman McInnes Notes: 432c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 433c96a6f78SLois Curfman McInnes 43436851e7fSLois Curfman McInnes Level: intermediate 43536851e7fSLois Curfman McInnes 4369b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4379b94acceSBarry Smith @*/ 43863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberUnsuccessfulSteps(SNES snes,PetscInt* nfails) 4399b94acceSBarry Smith { 4403a40ed3dSBarry Smith PetscFunctionBegin; 4414482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 4424482741eSBarry Smith PetscValidIntPointer(nfails,2); 44350ffb88aSMatthew Knepley *nfails = snes->numFailures; 44450ffb88aSMatthew Knepley PetscFunctionReturn(0); 44550ffb88aSMatthew Knepley } 44650ffb88aSMatthew Knepley 44750ffb88aSMatthew Knepley #undef __FUNCT__ 44850ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps" 44950ffb88aSMatthew Knepley /*@ 45050ffb88aSMatthew Knepley SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps 45150ffb88aSMatthew Knepley attempted by the nonlinear solver before it gives up. 45250ffb88aSMatthew Knepley 45350ffb88aSMatthew Knepley Not Collective 45450ffb88aSMatthew Knepley 45550ffb88aSMatthew Knepley Input Parameters: 45650ffb88aSMatthew Knepley + snes - SNES context 45750ffb88aSMatthew Knepley - maxFails - maximum of unsuccessful steps 45850ffb88aSMatthew Knepley 45950ffb88aSMatthew Knepley Level: intermediate 46050ffb88aSMatthew Knepley 46150ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 46250ffb88aSMatthew Knepley @*/ 46363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaximumUnsuccessfulSteps(SNES snes, PetscInt maxFails) 46450ffb88aSMatthew Knepley { 46550ffb88aSMatthew Knepley PetscFunctionBegin; 4664482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 46750ffb88aSMatthew Knepley snes->maxFailures = maxFails; 46850ffb88aSMatthew Knepley PetscFunctionReturn(0); 46950ffb88aSMatthew Knepley } 47050ffb88aSMatthew Knepley 47150ffb88aSMatthew Knepley #undef __FUNCT__ 47250ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps" 47350ffb88aSMatthew Knepley /*@ 47450ffb88aSMatthew Knepley SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps 47550ffb88aSMatthew Knepley attempted by the nonlinear solver before it gives up. 47650ffb88aSMatthew Knepley 47750ffb88aSMatthew Knepley Not Collective 47850ffb88aSMatthew Knepley 47950ffb88aSMatthew Knepley Input Parameter: 48050ffb88aSMatthew Knepley . snes - SNES context 48150ffb88aSMatthew Knepley 48250ffb88aSMatthew Knepley Output Parameter: 48350ffb88aSMatthew Knepley . maxFails - maximum of unsuccessful steps 48450ffb88aSMatthew Knepley 48550ffb88aSMatthew Knepley Level: intermediate 48650ffb88aSMatthew Knepley 48750ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 48850ffb88aSMatthew Knepley @*/ 48963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaximumUnsuccessfulSteps(SNES snes, PetscInt *maxFails) 49050ffb88aSMatthew Knepley { 49150ffb88aSMatthew Knepley PetscFunctionBegin; 4924482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 4934482741eSBarry Smith PetscValidIntPointer(maxFails,2); 49450ffb88aSMatthew Knepley *maxFails = snes->maxFailures; 4953a40ed3dSBarry Smith PetscFunctionReturn(0); 4969b94acceSBarry Smith } 497a847f771SSatish Balay 4984a2ae208SSatish Balay #undef __FUNCT__ 4994a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations" 500c96a6f78SLois Curfman McInnes /*@ 501c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 502c96a6f78SLois Curfman McInnes used by the nonlinear solver. 503c96a6f78SLois Curfman McInnes 504c7afd0dbSLois Curfman McInnes Not Collective 505c7afd0dbSLois Curfman McInnes 506c96a6f78SLois Curfman McInnes Input Parameter: 507c96a6f78SLois Curfman McInnes . snes - SNES context 508c96a6f78SLois Curfman McInnes 509c96a6f78SLois Curfman McInnes Output Parameter: 510c96a6f78SLois Curfman McInnes . lits - number of linear iterations 511c96a6f78SLois Curfman McInnes 512c96a6f78SLois Curfman McInnes Notes: 513c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 514c96a6f78SLois Curfman McInnes 51536851e7fSLois Curfman McInnes Level: intermediate 51636851e7fSLois Curfman McInnes 517c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 518*2b668275SBarry Smith 519*2b668275SBarry Smith .seealso: SNESGetIterationNumber(), SNESGetFunctionNorm() 520c96a6f78SLois Curfman McInnes @*/ 52163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberLinearIterations(SNES snes,PetscInt* lits) 522c96a6f78SLois Curfman McInnes { 5233a40ed3dSBarry Smith PetscFunctionBegin; 5244482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 5254482741eSBarry Smith PetscValidIntPointer(lits,2); 526c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 5273a40ed3dSBarry Smith PetscFunctionReturn(0); 528c96a6f78SLois Curfman McInnes } 529c96a6f78SLois Curfman McInnes 5304a2ae208SSatish Balay #undef __FUNCT__ 53194b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP" 53252baeb72SSatish Balay /*@ 53394b7f48cSBarry Smith SNESGetKSP - Returns the KSP context for a SNES solver. 5349b94acceSBarry Smith 53594b7f48cSBarry Smith Not Collective, but if SNES object is parallel, then KSP object is parallel 536c7afd0dbSLois Curfman McInnes 5379b94acceSBarry Smith Input Parameter: 5389b94acceSBarry Smith . snes - the SNES context 5399b94acceSBarry Smith 5409b94acceSBarry Smith Output Parameter: 54194b7f48cSBarry Smith . ksp - the KSP context 5429b94acceSBarry Smith 5439b94acceSBarry Smith Notes: 54494b7f48cSBarry Smith The user can then directly manipulate the KSP context to set various 5459b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5469b94acceSBarry Smith KSP and PC contexts as well. 5479b94acceSBarry Smith 54836851e7fSLois Curfman McInnes Level: beginner 54936851e7fSLois Curfman McInnes 55094b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context 5519b94acceSBarry Smith 55294b7f48cSBarry Smith .seealso: KSPGetPC() 5539b94acceSBarry Smith @*/ 55463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp) 5559b94acceSBarry Smith { 5563a40ed3dSBarry Smith PetscFunctionBegin; 5574482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 5584482741eSBarry Smith PetscValidPointer(ksp,2); 55994b7f48cSBarry Smith *ksp = snes->ksp; 5603a40ed3dSBarry Smith PetscFunctionReturn(0); 5619b94acceSBarry Smith } 56282bf6240SBarry Smith 5634a2ae208SSatish Balay #undef __FUNCT__ 5644a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc" 5656849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 566e24b481bSBarry Smith { 567e24b481bSBarry Smith PetscFunctionBegin; 568e24b481bSBarry Smith PetscFunctionReturn(0); 569e24b481bSBarry Smith } 570e24b481bSBarry Smith 5719b94acceSBarry Smith /* -----------------------------------------------------------*/ 5724a2ae208SSatish Balay #undef __FUNCT__ 5734a2ae208SSatish Balay #define __FUNCT__ "SNESCreate" 57452baeb72SSatish Balay /*@ 5759b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5769b94acceSBarry Smith 577c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 578c7afd0dbSLois Curfman McInnes 579c7afd0dbSLois Curfman McInnes Input Parameters: 580c7afd0dbSLois Curfman McInnes + comm - MPI communicator 5819b94acceSBarry Smith 5829b94acceSBarry Smith Output Parameter: 5839b94acceSBarry Smith . outsnes - the new SNES context 5849b94acceSBarry Smith 585c7afd0dbSLois Curfman McInnes Options Database Keys: 586c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 587c7afd0dbSLois Curfman McInnes and no preconditioning matrix 588c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 589c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 590c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 591c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 592c1f60f51SBarry Smith 59336851e7fSLois Curfman McInnes Level: beginner 59436851e7fSLois Curfman McInnes 5959b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5969b94acceSBarry Smith 5974b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES 5989b94acceSBarry Smith @*/ 59963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes) 6009b94acceSBarry Smith { 601dfbe8321SBarry Smith PetscErrorCode ierr; 6029b94acceSBarry Smith SNES snes; 6039b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 60437fcc0dbSBarry Smith 6053a40ed3dSBarry Smith PetscFunctionBegin; 606ed1caa07SMatthew Knepley PetscValidPointer(outsnes,2); 6078ba1e511SMatthew Knepley *outsnes = PETSC_NULL; 6088ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 6098ba1e511SMatthew Knepley ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 6108ba1e511SMatthew Knepley #endif 6118ba1e511SMatthew Knepley 61252e6d16bSBarry Smith ierr = PetscHeaderCreate(snes,_p_SNES,PetscInt,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 613e24b481bSBarry Smith snes->bops->publish = SNESPublish_Petsc; 6149b94acceSBarry Smith snes->max_its = 50; 6159750a799SBarry Smith snes->max_funcs = 10000; 6169b94acceSBarry Smith snes->norm = 0.0; 617b4874afaSBarry Smith snes->rtol = 1.e-8; 618b4874afaSBarry Smith snes->ttol = 0.0; 61970441072SBarry Smith snes->abstol = 1.e-50; 6209b94acceSBarry Smith snes->xtol = 1.e-8; 6214b27c08aSLois Curfman McInnes snes->deltatol = 1.e-12; 6229b94acceSBarry Smith snes->nfuncs = 0; 62350ffb88aSMatthew Knepley snes->numFailures = 0; 62450ffb88aSMatthew Knepley snes->maxFailures = 1; 6257a00f4a9SLois Curfman McInnes snes->linear_its = 0; 626639f9d9dSBarry Smith snes->numbermonitors = 0; 6279b94acceSBarry Smith snes->data = 0; 6289b94acceSBarry Smith snes->view = 0; 6294dc4c822SBarry Smith snes->setupcalled = PETSC_FALSE; 630186905e3SBarry Smith snes->ksp_ewconv = PETSC_FALSE; 6316f24a144SLois Curfman McInnes snes->vwork = 0; 6326f24a144SLois Curfman McInnes snes->nwork = 0; 633758f92a0SBarry Smith snes->conv_hist_len = 0; 634758f92a0SBarry Smith snes->conv_hist_max = 0; 635758f92a0SBarry Smith snes->conv_hist = PETSC_NULL; 636758f92a0SBarry Smith snes->conv_hist_its = PETSC_NULL; 637758f92a0SBarry Smith snes->conv_hist_reset = PETSC_TRUE; 638184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 6399b94acceSBarry Smith 6409b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 641b0a32e0cSBarry Smith ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 64252e6d16bSBarry Smith ierr = PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));CHKERRQ(ierr); 6439b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6449b94acceSBarry Smith kctx->version = 2; 6459b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6469b94acceSBarry Smith this was too large for some test cases */ 6479b94acceSBarry Smith kctx->rtol_last = 0; 6489b94acceSBarry Smith kctx->rtol_max = .9; 6499b94acceSBarry Smith kctx->gamma = 1.0; 6509b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6519b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6529b94acceSBarry Smith kctx->threshold = .1; 6539b94acceSBarry Smith kctx->lresid_last = 0; 6549b94acceSBarry Smith kctx->norm_last = 0; 6559b94acceSBarry Smith 65694b7f48cSBarry Smith ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr); 65752e6d16bSBarry Smith ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr); 6585334005bSBarry Smith 6599b94acceSBarry Smith *outsnes = snes; 66000036973SBarry Smith ierr = PetscPublishAll(snes);CHKERRQ(ierr); 6613a40ed3dSBarry Smith PetscFunctionReturn(0); 6629b94acceSBarry Smith } 6639b94acceSBarry Smith 6644a2ae208SSatish Balay #undef __FUNCT__ 6654a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction" 6669b94acceSBarry Smith /*@C 6679b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6689b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6699b94acceSBarry Smith equations. 6709b94acceSBarry Smith 671fee21e36SBarry Smith Collective on SNES 672fee21e36SBarry Smith 673c7afd0dbSLois Curfman McInnes Input Parameters: 674c7afd0dbSLois Curfman McInnes + snes - the SNES context 675c7afd0dbSLois Curfman McInnes . func - function evaluation routine 676c7afd0dbSLois Curfman McInnes . r - vector to store function value 677c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 678c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6799b94acceSBarry Smith 680c7afd0dbSLois Curfman McInnes Calling sequence of func: 6818d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 682c7afd0dbSLois Curfman McInnes 683313e4042SLois Curfman McInnes . f - function vector 684c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 6859b94acceSBarry Smith 6869b94acceSBarry Smith Notes: 6879b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6889b94acceSBarry Smith $ f'(x) x = -f(x), 689c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 6909b94acceSBarry Smith 69136851e7fSLois Curfman McInnes Level: beginner 69236851e7fSLois Curfman McInnes 6939b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6949b94acceSBarry Smith 695a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6969b94acceSBarry Smith @*/ 69763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 6989b94acceSBarry Smith { 6993a40ed3dSBarry Smith PetscFunctionBegin; 7004482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7014482741eSBarry Smith PetscValidHeaderSpecific(r,VEC_COOKIE,2); 702c9780b6fSBarry Smith PetscCheckSameComm(snes,1,r,2); 703184914b5SBarry Smith 7049b94acceSBarry Smith snes->computefunction = func; 7059b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7069b94acceSBarry Smith snes->funP = ctx; 7073a40ed3dSBarry Smith PetscFunctionReturn(0); 7089b94acceSBarry Smith } 7099b94acceSBarry Smith 7103ab0aad5SBarry Smith /* --------------------------------------------------------------- */ 7113ab0aad5SBarry Smith #undef __FUNCT__ 7123ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs" 7133ab0aad5SBarry Smith /*@C 7143ab0aad5SBarry Smith SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set 7153ab0aad5SBarry Smith it assumes a zero right hand side. 7163ab0aad5SBarry Smith 7173ab0aad5SBarry Smith Collective on SNES 7183ab0aad5SBarry Smith 7193ab0aad5SBarry Smith Input Parameters: 7203ab0aad5SBarry Smith + snes - the SNES context 7213ab0aad5SBarry Smith - rhs - the right hand side vector or PETSC_NULL for a zero right hand side 7223ab0aad5SBarry Smith 7233ab0aad5SBarry Smith Level: intermediate 7243ab0aad5SBarry Smith 7253ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side 7263ab0aad5SBarry Smith 7271096aae1SMatthew Knepley .seealso: SNESGetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 7283ab0aad5SBarry Smith @*/ 72963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetRhs(SNES snes,Vec rhs) 7303ab0aad5SBarry Smith { 731dfbe8321SBarry Smith PetscErrorCode ierr; 7323ab0aad5SBarry Smith 7333ab0aad5SBarry Smith PetscFunctionBegin; 7343ab0aad5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7353ab0aad5SBarry Smith if (rhs) { 7363ab0aad5SBarry Smith PetscValidHeaderSpecific(rhs,VEC_COOKIE,2); 7373ab0aad5SBarry Smith PetscCheckSameComm(snes,1,rhs,2); 7383ab0aad5SBarry Smith ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr); 7393ab0aad5SBarry Smith } 7403ab0aad5SBarry Smith if (snes->afine) { 7413ab0aad5SBarry Smith ierr = VecDestroy(snes->afine);CHKERRQ(ierr); 7423ab0aad5SBarry Smith } 7433ab0aad5SBarry Smith snes->afine = rhs; 7443ab0aad5SBarry Smith PetscFunctionReturn(0); 7453ab0aad5SBarry Smith } 7463ab0aad5SBarry Smith 7474a2ae208SSatish Balay #undef __FUNCT__ 7481096aae1SMatthew Knepley #define __FUNCT__ "SNESGetRhs" 7491096aae1SMatthew Knepley /*@C 7501096aae1SMatthew Knepley SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 7511096aae1SMatthew Knepley it assumes a zero right hand side. 7521096aae1SMatthew Knepley 7531096aae1SMatthew Knepley Collective on SNES 7541096aae1SMatthew Knepley 7551096aae1SMatthew Knepley Input Parameter: 7561096aae1SMatthew Knepley . snes - the SNES context 7571096aae1SMatthew Knepley 7581096aae1SMatthew Knepley Output Parameter: 7591096aae1SMatthew Knepley . rhs - the right hand side vector or PETSC_NULL for a zero right hand side 7601096aae1SMatthew Knepley 7611096aae1SMatthew Knepley Level: intermediate 7621096aae1SMatthew Knepley 7631096aae1SMatthew Knepley .keywords: SNES, nonlinear, get, function, right hand side 7641096aae1SMatthew Knepley 7651096aae1SMatthew Knepley .seealso: SNESSetRhs(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 7661096aae1SMatthew Knepley @*/ 7671096aae1SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs) 7681096aae1SMatthew Knepley { 7691096aae1SMatthew Knepley PetscFunctionBegin; 7701096aae1SMatthew Knepley PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7711096aae1SMatthew Knepley PetscValidPointer(rhs,2); 7721096aae1SMatthew Knepley *rhs = snes->afine; 7731096aae1SMatthew Knepley PetscFunctionReturn(0); 7741096aae1SMatthew Knepley } 7751096aae1SMatthew Knepley 7761096aae1SMatthew Knepley #undef __FUNCT__ 7774a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction" 7789b94acceSBarry Smith /*@ 77936851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7809b94acceSBarry Smith SNESSetFunction(). 7819b94acceSBarry Smith 782c7afd0dbSLois Curfman McInnes Collective on SNES 783c7afd0dbSLois Curfman McInnes 7849b94acceSBarry Smith Input Parameters: 785c7afd0dbSLois Curfman McInnes + snes - the SNES context 786c7afd0dbSLois Curfman McInnes - x - input vector 7879b94acceSBarry Smith 7889b94acceSBarry Smith Output Parameter: 7893638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7909b94acceSBarry Smith 7911bffabb2SLois Curfman McInnes Notes: 79236851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 79336851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 79436851e7fSLois Curfman McInnes themselves. 79536851e7fSLois Curfman McInnes 79636851e7fSLois Curfman McInnes Level: developer 79736851e7fSLois Curfman McInnes 7989b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 7999b94acceSBarry Smith 800a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 8019b94acceSBarry Smith @*/ 80263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y) 8039b94acceSBarry Smith { 804dfbe8321SBarry Smith PetscErrorCode ierr; 8059b94acceSBarry Smith 8063a40ed3dSBarry Smith PetscFunctionBegin; 8074482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8084482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 8094482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,3); 810c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 811c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,3); 812184914b5SBarry Smith 813d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 8141096aae1SMatthew Knepley if (snes->computefunction) { 815d64ed03dSBarry Smith PetscStackPush("SNES user function"); 816e9a2bbcdSBarry Smith CHKMEMQ; 81719717074SBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); 818e9a2bbcdSBarry Smith CHKMEMQ; 819d64ed03dSBarry Smith PetscStackPop; 820d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 82119717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 82219717074SBarry Smith } 823d5e45103SBarry Smith CHKERRQ(ierr); 8241096aae1SMatthew Knepley } else if (snes->afine) { 8251096aae1SMatthew Knepley ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 8261096aae1SMatthew Knepley } else { 8271096aae1SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve()."); 8281096aae1SMatthew Knepley } 8293ab0aad5SBarry Smith if (snes->afine) { 830016dedfbSBarry Smith ierr = VecAXPY(y,-1.0,snes->afine);CHKERRQ(ierr); 8313ab0aad5SBarry Smith } 832ae3c334cSLois Curfman McInnes snes->nfuncs++; 833d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 8343a40ed3dSBarry Smith PetscFunctionReturn(0); 8359b94acceSBarry Smith } 8369b94acceSBarry Smith 8374a2ae208SSatish Balay #undef __FUNCT__ 8384a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian" 83962fef451SLois Curfman McInnes /*@ 84062fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 84162fef451SLois Curfman McInnes set with SNESSetJacobian(). 84262fef451SLois Curfman McInnes 843c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 844c7afd0dbSLois Curfman McInnes 84562fef451SLois Curfman McInnes Input Parameters: 846c7afd0dbSLois Curfman McInnes + snes - the SNES context 847c7afd0dbSLois Curfman McInnes - x - input vector 84862fef451SLois Curfman McInnes 84962fef451SLois Curfman McInnes Output Parameters: 850c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 85162fef451SLois Curfman McInnes . B - optional preconditioning matrix 852*2b668275SBarry Smith - flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER) 853fee21e36SBarry Smith 85462fef451SLois Curfman McInnes Notes: 85562fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 85662fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 85762fef451SLois Curfman McInnes 85894b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 859dc5a77f8SLois Curfman McInnes flag parameter. 86062fef451SLois Curfman McInnes 86136851e7fSLois Curfman McInnes Level: developer 86236851e7fSLois Curfman McInnes 86362fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 86462fef451SLois Curfman McInnes 865*2b668275SBarry Smith .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure 86662fef451SLois Curfman McInnes @*/ 86763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8689b94acceSBarry Smith { 869dfbe8321SBarry Smith PetscErrorCode ierr; 8703a40ed3dSBarry Smith 8713a40ed3dSBarry Smith PetscFunctionBegin; 8724482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8734482741eSBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE,2); 8744482741eSBarry Smith PetscValidPointer(flg,5); 875c9780b6fSBarry Smith PetscCheckSameComm(snes,1,X,2); 8763a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 877d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 878c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 879d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 880dc67937bSBarry Smith CHKMEMQ; 8819b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 882dc67937bSBarry Smith CHKMEMQ; 883d64ed03dSBarry Smith PetscStackPop; 884d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 8856d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 8864482741eSBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 8874482741eSBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 8883a40ed3dSBarry Smith PetscFunctionReturn(0); 8899b94acceSBarry Smith } 8909b94acceSBarry Smith 8914a2ae208SSatish Balay #undef __FUNCT__ 8924a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian" 8939b94acceSBarry Smith /*@C 8949b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 895044dda88SLois Curfman McInnes location to store the matrix. 8969b94acceSBarry Smith 897c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 898c7afd0dbSLois Curfman McInnes 8999b94acceSBarry Smith Input Parameters: 900c7afd0dbSLois Curfman McInnes + snes - the SNES context 9019b94acceSBarry Smith . A - Jacobian matrix 9029b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 9039b94acceSBarry Smith . func - Jacobian evaluation routine 904c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 9052cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 9069b94acceSBarry Smith 9079b94acceSBarry Smith Calling sequence of func: 9088d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9099b94acceSBarry Smith 910c7afd0dbSLois Curfman McInnes + x - input vector 9119b94acceSBarry Smith . A - Jacobian matrix 9129b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 913ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 914*2b668275SBarry Smith structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER 915c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 9169b94acceSBarry Smith 9179b94acceSBarry Smith Notes: 91894b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 9192cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 920ac21db08SLois Curfman McInnes 921ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9229b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9239b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9249b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9259b94acceSBarry Smith throughout the global iterations. 9269b94acceSBarry Smith 92736851e7fSLois Curfman McInnes Level: beginner 92836851e7fSLois Curfman McInnes 9299b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9309b94acceSBarry Smith 931*2b668275SBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure 9329b94acceSBarry Smith @*/ 93363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 9349b94acceSBarry Smith { 935dfbe8321SBarry Smith PetscErrorCode ierr; 9363a7fca6bSBarry Smith 9373a40ed3dSBarry Smith PetscFunctionBegin; 9384482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 9394482741eSBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 9404482741eSBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 941c9780b6fSBarry Smith if (A) PetscCheckSameComm(snes,1,A,2); 942c9780b6fSBarry Smith if (B) PetscCheckSameComm(snes,1,B,2); 9433a7fca6bSBarry Smith if (func) snes->computejacobian = func; 9443a7fca6bSBarry Smith if (ctx) snes->jacP = ctx; 9453a7fca6bSBarry Smith if (A) { 9463a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 9479b94acceSBarry Smith snes->jacobian = A; 9483a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 9493a7fca6bSBarry Smith } 9503a7fca6bSBarry Smith if (B) { 9513a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 9529b94acceSBarry Smith snes->jacobian_pre = B; 9533a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 9543a7fca6bSBarry Smith } 95569a612a9SBarry Smith ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 9563a40ed3dSBarry Smith PetscFunctionReturn(0); 9579b94acceSBarry Smith } 95862fef451SLois Curfman McInnes 9594a2ae208SSatish Balay #undef __FUNCT__ 9604a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian" 961c2aafc4cSSatish Balay /*@C 962b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 963b4fd4287SBarry Smith provided context for evaluating the Jacobian. 964b4fd4287SBarry Smith 965c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 966c7afd0dbSLois Curfman McInnes 967b4fd4287SBarry Smith Input Parameter: 968b4fd4287SBarry Smith . snes - the nonlinear solver context 969b4fd4287SBarry Smith 970b4fd4287SBarry Smith Output Parameters: 971c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 972b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 97370e92668SMatthew Knepley . func - location to put Jacobian function (or PETSC_NULL) 97470e92668SMatthew Knepley - ctx - location to stash Jacobian ctx (or PETSC_NULL) 975fee21e36SBarry Smith 97636851e7fSLois Curfman McInnes Level: advanced 97736851e7fSLois Curfman McInnes 978b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 979b4fd4287SBarry Smith @*/ 98070e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 981b4fd4287SBarry Smith { 9823a40ed3dSBarry Smith PetscFunctionBegin; 9834482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 984b4fd4287SBarry Smith if (A) *A = snes->jacobian; 985b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 98600036973SBarry Smith if (func) *func = snes->computejacobian; 98770e92668SMatthew Knepley if (ctx) *ctx = snes->jacP; 9883a40ed3dSBarry Smith PetscFunctionReturn(0); 989b4fd4287SBarry Smith } 990b4fd4287SBarry Smith 9919b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 99263dd3a1aSKris Buschelman EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 9939b94acceSBarry Smith 9944a2ae208SSatish Balay #undef __FUNCT__ 9954a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp" 9969b94acceSBarry Smith /*@ 9979b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 998272ac6f2SLois Curfman McInnes of a nonlinear solver. 9999b94acceSBarry Smith 1000fee21e36SBarry Smith Collective on SNES 1001fee21e36SBarry Smith 1002c7afd0dbSLois Curfman McInnes Input Parameters: 100370e92668SMatthew Knepley . snes - the SNES context 1004c7afd0dbSLois Curfman McInnes 1005272ac6f2SLois Curfman McInnes Notes: 1006272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1007272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1008272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1009272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1010272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1011272ac6f2SLois Curfman McInnes 101236851e7fSLois Curfman McInnes Level: advanced 101336851e7fSLois Curfman McInnes 10149b94acceSBarry Smith .keywords: SNES, nonlinear, setup 10159b94acceSBarry Smith 10169b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 10179b94acceSBarry Smith @*/ 101870e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 10199b94acceSBarry Smith { 1020dfbe8321SBarry Smith PetscErrorCode ierr; 10214b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 10223a40ed3dSBarry Smith 10233a40ed3dSBarry Smith PetscFunctionBegin; 10244482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 10254dc4c822SBarry Smith if (snes->setupcalled) PetscFunctionReturn(0); 10269b94acceSBarry Smith 1027b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1028c1f60f51SBarry Smith /* 1029c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1030dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1031c1f60f51SBarry Smith */ 1032c1f60f51SBarry Smith if (flg) { 1033c1f60f51SBarry Smith Mat J; 10345a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10355a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1036ae15b995SBarry Smith ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr); 10373a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 10383a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 1039c1f60f51SBarry Smith } 104045fc7adcSBarry Smith 104103c60df9SBarry Smith #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT) 104245fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 104345fc7adcSBarry Smith if (flg) { 104445fc7adcSBarry Smith Mat J; 104545fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 104645fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 104745fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 104845fc7adcSBarry Smith } 104932a4b47aSMatthew Knepley #endif 105045fc7adcSBarry Smith 1051b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1052c1f60f51SBarry Smith /* 1053dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1054c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1055c1f60f51SBarry Smith */ 105631615d04SLois Curfman McInnes if (flg) { 1057272ac6f2SLois Curfman McInnes Mat J; 1058b5d62d44SBarry Smith KSP ksp; 105994b7f48cSBarry Smith PC pc; 1060f3ef73ceSBarry Smith 10615a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10623a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1063ae15b995SBarry Smith ierr = PetscInfo(snes,"Setting default matrix-free operator and preconditioner routines;\nThat is no preconditioner is being used.\n");CHKERRQ(ierr); 10643a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 10653a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 10663a7fca6bSBarry Smith 1067f3ef73ceSBarry Smith /* force no preconditioner */ 106894b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1069b5d62d44SBarry Smith ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1070a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1071a9815358SBarry Smith if (!flg) { 1072f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1073272ac6f2SLois Curfman McInnes } 1074a9815358SBarry Smith } 1075f3ef73ceSBarry Smith 10761096aae1SMatthew Knepley if (!snes->vec_func && !snes->afine) { 10771096aae1SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 10781096aae1SMatthew Knepley } 10791096aae1SMatthew Knepley if (!snes->computefunction && !snes->afine) { 10801096aae1SMatthew Knepley SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 10811096aae1SMatthew Knepley } 108229bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1083a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 108429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1085a8c6a408SBarry Smith } 1086a703fe33SLois Curfman McInnes 1087a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 10884b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 10896831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 109094b7f48cSBarry Smith KSP ksp; 109194b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1092186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1093a703fe33SLois Curfman McInnes } 10944b27c08aSLois Curfman McInnes 1095a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 10967aaed0d8SBarry Smith snes->setupcalled = PETSC_TRUE; 10973a40ed3dSBarry Smith PetscFunctionReturn(0); 10989b94acceSBarry Smith } 10999b94acceSBarry Smith 11004a2ae208SSatish Balay #undef __FUNCT__ 11014a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 110252baeb72SSatish Balay /*@ 11039b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 11049b94acceSBarry Smith with SNESCreate(). 11059b94acceSBarry Smith 1106c7afd0dbSLois Curfman McInnes Collective on SNES 1107c7afd0dbSLois Curfman McInnes 11089b94acceSBarry Smith Input Parameter: 11099b94acceSBarry Smith . snes - the SNES context 11109b94acceSBarry Smith 111136851e7fSLois Curfman McInnes Level: beginner 111236851e7fSLois Curfman McInnes 11139b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 11149b94acceSBarry Smith 111563a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 11169b94acceSBarry Smith @*/ 111763dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 11189b94acceSBarry Smith { 11196849ba73SBarry Smith PetscErrorCode ierr; 11203a40ed3dSBarry Smith 11213a40ed3dSBarry Smith PetscFunctionBegin; 11224482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 11233a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1124d4bb536fSBarry Smith 1125be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 11260f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1127be0abb6dSBarry Smith 1128e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 112905b42c5fSBarry Smith ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr); 11303a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 11313a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 11323ab0aad5SBarry Smith if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 113394b7f48cSBarry Smith ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1134522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1135d952e501SBarry Smith ierr = SNESClearMonitor(snes);CHKERRQ(ierr); 1136a79aaaedSSatish Balay ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 11373a40ed3dSBarry Smith PetscFunctionReturn(0); 11389b94acceSBarry Smith } 11399b94acceSBarry Smith 11409b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 11419b94acceSBarry Smith 11424a2ae208SSatish Balay #undef __FUNCT__ 11434a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 11449b94acceSBarry Smith /*@ 1145d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 11469b94acceSBarry Smith 1147c7afd0dbSLois Curfman McInnes Collective on SNES 1148c7afd0dbSLois Curfman McInnes 11499b94acceSBarry Smith Input Parameters: 1150c7afd0dbSLois Curfman McInnes + snes - the SNES context 115170441072SBarry Smith . abstol - absolute convergence tolerance 115233174efeSLois Curfman McInnes . rtol - relative convergence tolerance 115333174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 115433174efeSLois Curfman McInnes of the change in the solution between steps 115533174efeSLois Curfman McInnes . maxit - maximum number of iterations 1156c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1157fee21e36SBarry Smith 115833174efeSLois Curfman McInnes Options Database Keys: 115970441072SBarry Smith + -snes_atol <abstol> - Sets abstol 1160c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1161c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1162c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1163c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 11649b94acceSBarry Smith 1165d7a720efSLois Curfman McInnes Notes: 11669b94acceSBarry Smith The default maximum number of iterations is 50. 11679b94acceSBarry Smith The default maximum number of function evaluations is 1000. 11689b94acceSBarry Smith 116936851e7fSLois Curfman McInnes Level: intermediate 117036851e7fSLois Curfman McInnes 117133174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 11729b94acceSBarry Smith 11732492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance() 11749b94acceSBarry Smith @*/ 117563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 11769b94acceSBarry Smith { 11773a40ed3dSBarry Smith PetscFunctionBegin; 11784482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 117970441072SBarry Smith if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1180d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1181d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1182d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1183d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11843a40ed3dSBarry Smith PetscFunctionReturn(0); 11859b94acceSBarry Smith } 11869b94acceSBarry Smith 11874a2ae208SSatish Balay #undef __FUNCT__ 11884a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11899b94acceSBarry Smith /*@ 119033174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 119133174efeSLois Curfman McInnes 1192c7afd0dbSLois Curfman McInnes Not Collective 1193c7afd0dbSLois Curfman McInnes 119433174efeSLois Curfman McInnes Input Parameters: 1195c7afd0dbSLois Curfman McInnes + snes - the SNES context 119670441072SBarry Smith . abstol - absolute convergence tolerance 119733174efeSLois Curfman McInnes . rtol - relative convergence tolerance 119833174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 119933174efeSLois Curfman McInnes of the change in the solution between steps 120033174efeSLois Curfman McInnes . maxit - maximum number of iterations 1201c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1202fee21e36SBarry Smith 120333174efeSLois Curfman McInnes Notes: 120433174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 120533174efeSLois Curfman McInnes 120636851e7fSLois Curfman McInnes Level: intermediate 120736851e7fSLois Curfman McInnes 120833174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 120933174efeSLois Curfman McInnes 121033174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 121133174efeSLois Curfman McInnes @*/ 121263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 121333174efeSLois Curfman McInnes { 12143a40ed3dSBarry Smith PetscFunctionBegin; 12154482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 121670441072SBarry Smith if (abstol) *abstol = snes->abstol; 121733174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 121833174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 121933174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 122033174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 12213a40ed3dSBarry Smith PetscFunctionReturn(0); 122233174efeSLois Curfman McInnes } 122333174efeSLois Curfman McInnes 12244a2ae208SSatish Balay #undef __FUNCT__ 12254a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 122633174efeSLois Curfman McInnes /*@ 12279b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 12289b94acceSBarry Smith 1229fee21e36SBarry Smith Collective on SNES 1230fee21e36SBarry Smith 1231c7afd0dbSLois Curfman McInnes Input Parameters: 1232c7afd0dbSLois Curfman McInnes + snes - the SNES context 1233c7afd0dbSLois Curfman McInnes - tol - tolerance 1234c7afd0dbSLois Curfman McInnes 12359b94acceSBarry Smith Options Database Key: 1236c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 12379b94acceSBarry Smith 123836851e7fSLois Curfman McInnes Level: intermediate 123936851e7fSLois Curfman McInnes 12409b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 12419b94acceSBarry Smith 12422492ecdbSBarry Smith .seealso: SNESSetTolerances() 12439b94acceSBarry Smith @*/ 124463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 12459b94acceSBarry Smith { 12463a40ed3dSBarry Smith PetscFunctionBegin; 12474482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 12489b94acceSBarry Smith snes->deltatol = tol; 12493a40ed3dSBarry Smith PetscFunctionReturn(0); 12509b94acceSBarry Smith } 12519b94acceSBarry Smith 1252df9fa365SBarry Smith /* 1253df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1254df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1255df9fa365SBarry Smith macros instead of functions 1256df9fa365SBarry Smith */ 12574a2ae208SSatish Balay #undef __FUNCT__ 12584a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 125963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1260ce1608b8SBarry Smith { 1261dfbe8321SBarry Smith PetscErrorCode ierr; 1262ce1608b8SBarry Smith 1263ce1608b8SBarry Smith PetscFunctionBegin; 12644482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1265ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1266ce1608b8SBarry Smith PetscFunctionReturn(0); 1267ce1608b8SBarry Smith } 1268ce1608b8SBarry Smith 12694a2ae208SSatish Balay #undef __FUNCT__ 12704a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 127163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1272df9fa365SBarry Smith { 1273dfbe8321SBarry Smith PetscErrorCode ierr; 1274df9fa365SBarry Smith 1275df9fa365SBarry Smith PetscFunctionBegin; 1276df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1277df9fa365SBarry Smith PetscFunctionReturn(0); 1278df9fa365SBarry Smith } 1279df9fa365SBarry Smith 12804a2ae208SSatish Balay #undef __FUNCT__ 12814a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 128263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESLGMonitorDestroy(PetscDrawLG draw) 1283df9fa365SBarry Smith { 1284dfbe8321SBarry Smith PetscErrorCode ierr; 1285df9fa365SBarry Smith 1286df9fa365SBarry Smith PetscFunctionBegin; 1287df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1288df9fa365SBarry Smith PetscFunctionReturn(0); 1289df9fa365SBarry Smith } 1290df9fa365SBarry Smith 12919b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12929b94acceSBarry Smith 12934a2ae208SSatish Balay #undef __FUNCT__ 12944a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12959b94acceSBarry Smith /*@C 12965cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12979b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12989b94acceSBarry Smith progress. 12999b94acceSBarry Smith 1300fee21e36SBarry Smith Collective on SNES 1301fee21e36SBarry Smith 1302c7afd0dbSLois Curfman McInnes Input Parameters: 1303c7afd0dbSLois Curfman McInnes + snes - the SNES context 1304c7afd0dbSLois Curfman McInnes . func - monitoring routine 1305b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1306e8105e01SRichard Katz monitor routine (use PETSC_NULL if no context is desired) 1307b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1308b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 13099b94acceSBarry Smith 1310c7afd0dbSLois Curfman McInnes Calling sequence of func: 1311a7cc72afSBarry Smith $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1312c7afd0dbSLois Curfman McInnes 1313c7afd0dbSLois Curfman McInnes + snes - the SNES context 1314c7afd0dbSLois Curfman McInnes . its - iteration number 1315c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 131640a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 13179b94acceSBarry Smith 13189665c990SLois Curfman McInnes Options Database Keys: 1319c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1320c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1321c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1322c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1323c7afd0dbSLois Curfman McInnes been hardwired into a code by 1324c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1325c7afd0dbSLois Curfman McInnes does not cancel those set via 1326c7afd0dbSLois Curfman McInnes the options database. 13279665c990SLois Curfman McInnes 1328639f9d9dSBarry Smith Notes: 13296bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 13306bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 13316bc08f3fSLois Curfman McInnes order in which they were set. 1332639f9d9dSBarry Smith 133336851e7fSLois Curfman McInnes Level: intermediate 133436851e7fSLois Curfman McInnes 13359b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 13369b94acceSBarry Smith 13375cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 13389b94acceSBarry Smith @*/ 133963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 13409b94acceSBarry Smith { 13413a40ed3dSBarry Smith PetscFunctionBegin; 13424482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1343639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 134429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1345639f9d9dSBarry Smith } 1346639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1347b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1348639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 13493a40ed3dSBarry Smith PetscFunctionReturn(0); 13509b94acceSBarry Smith } 13519b94acceSBarry Smith 13524a2ae208SSatish Balay #undef __FUNCT__ 13534a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 13545cd90555SBarry Smith /*@C 13555cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 13565cd90555SBarry Smith 1357c7afd0dbSLois Curfman McInnes Collective on SNES 1358c7afd0dbSLois Curfman McInnes 13595cd90555SBarry Smith Input Parameters: 13605cd90555SBarry Smith . snes - the SNES context 13615cd90555SBarry Smith 13621a480d89SAdministrator Options Database Key: 1363c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1364c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1365c7afd0dbSLois Curfman McInnes set via the options database 13665cd90555SBarry Smith 13675cd90555SBarry Smith Notes: 13685cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 13695cd90555SBarry Smith 137036851e7fSLois Curfman McInnes Level: intermediate 137136851e7fSLois Curfman McInnes 13725cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 13735cd90555SBarry Smith 13745cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 13755cd90555SBarry Smith @*/ 137663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESClearMonitor(SNES snes) 13775cd90555SBarry Smith { 1378d952e501SBarry Smith PetscErrorCode ierr; 1379d952e501SBarry Smith PetscInt i; 1380d952e501SBarry Smith 13815cd90555SBarry Smith PetscFunctionBegin; 13824482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1383d952e501SBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1384d952e501SBarry Smith if (snes->monitordestroy[i]) { 1385d952e501SBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1386d952e501SBarry Smith } 1387d952e501SBarry Smith } 13885cd90555SBarry Smith snes->numbermonitors = 0; 13895cd90555SBarry Smith PetscFunctionReturn(0); 13905cd90555SBarry Smith } 13915cd90555SBarry Smith 13924a2ae208SSatish Balay #undef __FUNCT__ 13934a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13949b94acceSBarry Smith /*@C 13959b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13969b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13979b94acceSBarry Smith 1398fee21e36SBarry Smith Collective on SNES 1399fee21e36SBarry Smith 1400c7afd0dbSLois Curfman McInnes Input Parameters: 1401c7afd0dbSLois Curfman McInnes + snes - the SNES context 1402c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1403c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1404c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 14059b94acceSBarry Smith 1406c7afd0dbSLois Curfman McInnes Calling sequence of func: 14076849ba73SBarry Smith $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1408c7afd0dbSLois Curfman McInnes 1409c7afd0dbSLois Curfman McInnes + snes - the SNES context 1410c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1411184914b5SBarry Smith . reason - reason for convergence/divergence 1412c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 14134b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 14144b27c08aSLois Curfman McInnes - f - 2-norm of function 14159b94acceSBarry Smith 141636851e7fSLois Curfman McInnes Level: advanced 141736851e7fSLois Curfman McInnes 14189b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 14199b94acceSBarry Smith 14204b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 14219b94acceSBarry Smith @*/ 142263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 14239b94acceSBarry Smith { 14243a40ed3dSBarry Smith PetscFunctionBegin; 14254482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14269b94acceSBarry Smith (snes)->converged = func; 14279b94acceSBarry Smith (snes)->cnvP = cctx; 14283a40ed3dSBarry Smith PetscFunctionReturn(0); 14299b94acceSBarry Smith } 14309b94acceSBarry Smith 14314a2ae208SSatish Balay #undef __FUNCT__ 14324a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 143352baeb72SSatish Balay /*@ 1434184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1435184914b5SBarry Smith 1436184914b5SBarry Smith Not Collective 1437184914b5SBarry Smith 1438184914b5SBarry Smith Input Parameter: 1439184914b5SBarry Smith . snes - the SNES context 1440184914b5SBarry Smith 1441184914b5SBarry Smith Output Parameter: 1442e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1443184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1444184914b5SBarry Smith 1445184914b5SBarry Smith Level: intermediate 1446184914b5SBarry Smith 1447184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1448184914b5SBarry Smith 1449184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1450184914b5SBarry Smith 14514b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1452184914b5SBarry Smith @*/ 145363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1454184914b5SBarry Smith { 1455184914b5SBarry Smith PetscFunctionBegin; 14564482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14574482741eSBarry Smith PetscValidPointer(reason,2); 1458184914b5SBarry Smith *reason = snes->reason; 1459184914b5SBarry Smith PetscFunctionReturn(0); 1460184914b5SBarry Smith } 1461184914b5SBarry Smith 14624a2ae208SSatish Balay #undef __FUNCT__ 14634a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1464c9005455SLois Curfman McInnes /*@ 1465c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1466c9005455SLois Curfman McInnes 1467fee21e36SBarry Smith Collective on SNES 1468fee21e36SBarry Smith 1469c7afd0dbSLois Curfman McInnes Input Parameters: 1470c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1471c7afd0dbSLois Curfman McInnes . a - array to hold history 1472cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1473758f92a0SBarry Smith . na - size of a and its 147464731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1475758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1476c7afd0dbSLois Curfman McInnes 1477c9005455SLois Curfman McInnes Notes: 14784b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1479c9005455SLois Curfman McInnes at each step. 1480c9005455SLois Curfman McInnes 1481c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1482c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1483c9005455SLois Curfman McInnes during the section of code that is being timed. 1484c9005455SLois Curfman McInnes 148536851e7fSLois Curfman McInnes Level: intermediate 148636851e7fSLois Curfman McInnes 1487c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1488758f92a0SBarry Smith 148908405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1490758f92a0SBarry Smith 1491c9005455SLois Curfman McInnes @*/ 149263dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1493c9005455SLois Curfman McInnes { 14943a40ed3dSBarry Smith PetscFunctionBegin; 14954482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14964482741eSBarry Smith if (na) PetscValidScalarPointer(a,2); 1497c9005455SLois Curfman McInnes snes->conv_hist = a; 1498758f92a0SBarry Smith snes->conv_hist_its = its; 1499758f92a0SBarry Smith snes->conv_hist_max = na; 1500758f92a0SBarry Smith snes->conv_hist_reset = reset; 1501758f92a0SBarry Smith PetscFunctionReturn(0); 1502758f92a0SBarry Smith } 1503758f92a0SBarry Smith 15044a2ae208SSatish Balay #undef __FUNCT__ 15054a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 15060c4c9dddSBarry Smith /*@C 1507758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1508758f92a0SBarry Smith 1509758f92a0SBarry Smith Collective on SNES 1510758f92a0SBarry Smith 1511758f92a0SBarry Smith Input Parameter: 1512758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1513758f92a0SBarry Smith 1514758f92a0SBarry Smith Output Parameters: 1515758f92a0SBarry Smith . a - array to hold history 1516758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1517758f92a0SBarry Smith negative if not converged) for each solve. 1518758f92a0SBarry Smith - na - size of a and its 1519758f92a0SBarry Smith 1520758f92a0SBarry Smith Notes: 1521758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1522758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1523758f92a0SBarry Smith 1524758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1525758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1526758f92a0SBarry Smith during the section of code that is being timed. 1527758f92a0SBarry Smith 1528758f92a0SBarry Smith Level: intermediate 1529758f92a0SBarry Smith 1530758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1531758f92a0SBarry Smith 1532758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1533758f92a0SBarry Smith 1534758f92a0SBarry Smith @*/ 153563dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1536758f92a0SBarry Smith { 1537758f92a0SBarry Smith PetscFunctionBegin; 15384482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1539758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1540758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1541758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 15423a40ed3dSBarry Smith PetscFunctionReturn(0); 1543c9005455SLois Curfman McInnes } 1544c9005455SLois Curfman McInnes 1545e74ef692SMatthew Knepley #undef __FUNCT__ 1546e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 1547ac226902SBarry Smith /*@C 154876b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 15497e4bb74cSBarry Smith at the beginning o every iteration of the nonlinear solve. Specifically 15507e4bb74cSBarry Smith it is called just before the Jacobian is "evaluated". 155176b2cf59SMatthew Knepley 155276b2cf59SMatthew Knepley Collective on SNES 155376b2cf59SMatthew Knepley 155476b2cf59SMatthew Knepley Input Parameters: 155576b2cf59SMatthew Knepley . snes - The nonlinear solver context 155676b2cf59SMatthew Knepley . func - The function 155776b2cf59SMatthew Knepley 155876b2cf59SMatthew Knepley Calling sequence of func: 1559b5d30489SBarry Smith . func (SNES snes, PetscInt step); 156076b2cf59SMatthew Knepley 156176b2cf59SMatthew Knepley . step - The current step of the iteration 156276b2cf59SMatthew Knepley 156376b2cf59SMatthew Knepley Level: intermediate 156476b2cf59SMatthew Knepley 156576b2cf59SMatthew Knepley .keywords: SNES, update 1566b5d30489SBarry Smith 156776b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 156876b2cf59SMatthew Knepley @*/ 156963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 157076b2cf59SMatthew Knepley { 157176b2cf59SMatthew Knepley PetscFunctionBegin; 15724482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 157376b2cf59SMatthew Knepley snes->update = func; 157476b2cf59SMatthew Knepley PetscFunctionReturn(0); 157576b2cf59SMatthew Knepley } 157676b2cf59SMatthew Knepley 1577e74ef692SMatthew Knepley #undef __FUNCT__ 1578e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 157976b2cf59SMatthew Knepley /*@ 158076b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 158176b2cf59SMatthew Knepley 158276b2cf59SMatthew Knepley Not collective 158376b2cf59SMatthew Knepley 158476b2cf59SMatthew Knepley Input Parameters: 158576b2cf59SMatthew Knepley . snes - The nonlinear solver context 158676b2cf59SMatthew Knepley . step - The current step of the iteration 158776b2cf59SMatthew Knepley 1588205452f4SMatthew Knepley Level: intermediate 1589205452f4SMatthew Knepley 159076b2cf59SMatthew Knepley .keywords: SNES, update 159176b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 159276b2cf59SMatthew Knepley @*/ 159363dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 159476b2cf59SMatthew Knepley { 159576b2cf59SMatthew Knepley PetscFunctionBegin; 159676b2cf59SMatthew Knepley PetscFunctionReturn(0); 159776b2cf59SMatthew Knepley } 159876b2cf59SMatthew Knepley 15994a2ae208SSatish Balay #undef __FUNCT__ 16004a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 16019b94acceSBarry Smith /* 16029b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 16039b94acceSBarry Smith positive parameter delta. 16049b94acceSBarry Smith 16059b94acceSBarry Smith Input Parameters: 1606c7afd0dbSLois Curfman McInnes + snes - the SNES context 16079b94acceSBarry Smith . y - approximate solution of linear system 16089b94acceSBarry Smith . fnorm - 2-norm of current function 1609c7afd0dbSLois Curfman McInnes - delta - trust region size 16109b94acceSBarry Smith 16119b94acceSBarry Smith Output Parameters: 1612c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 16139b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 16149b94acceSBarry Smith region, and exceeds zero otherwise. 1615c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 16169b94acceSBarry Smith 16179b94acceSBarry Smith Note: 16184b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 16199b94acceSBarry Smith is set to be the maximum allowable step size. 16209b94acceSBarry Smith 16219b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16229b94acceSBarry Smith */ 1623dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 16249b94acceSBarry Smith { 1625064f8208SBarry Smith PetscReal nrm; 1626ea709b57SSatish Balay PetscScalar cnorm; 1627dfbe8321SBarry Smith PetscErrorCode ierr; 16283a40ed3dSBarry Smith 16293a40ed3dSBarry Smith PetscFunctionBegin; 16304482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16314482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1632c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,2); 1633184914b5SBarry Smith 1634064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1635064f8208SBarry Smith if (nrm > *delta) { 1636064f8208SBarry Smith nrm = *delta/nrm; 1637064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1638064f8208SBarry Smith cnorm = nrm; 16392dcb1b2aSMatthew Knepley ierr = VecScale(y,cnorm);CHKERRQ(ierr); 16409b94acceSBarry Smith *ynorm = *delta; 16419b94acceSBarry Smith } else { 16429b94acceSBarry Smith *gpnorm = 0.0; 1643064f8208SBarry Smith *ynorm = nrm; 16449b94acceSBarry Smith } 16453a40ed3dSBarry Smith PetscFunctionReturn(0); 16469b94acceSBarry Smith } 16479b94acceSBarry Smith 16484a2ae208SSatish Balay #undef __FUNCT__ 16494a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 16509b94acceSBarry Smith /*@ 1651f69a0ea3SMatthew Knepley SNESSolve - Solves a nonlinear system F(x) = b. 1652f69a0ea3SMatthew Knepley Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 16539b94acceSBarry Smith 1654c7afd0dbSLois Curfman McInnes Collective on SNES 1655c7afd0dbSLois Curfman McInnes 1656b2002411SLois Curfman McInnes Input Parameters: 1657c7afd0dbSLois Curfman McInnes + snes - the SNES context 1658f69a0ea3SMatthew Knepley . b - the constant part of the equation, or PETSC_NULL to use zero. 165970e92668SMatthew Knepley - x - the solution vector, or PETSC_NULL if it was set with SNESSetSolution() 16609b94acceSBarry Smith 1661b2002411SLois Curfman McInnes Notes: 16628ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 16638ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 16648ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 16658ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 16668ddd3da0SLois Curfman McInnes 166736851e7fSLois Curfman McInnes Level: beginner 166836851e7fSLois Curfman McInnes 16699b94acceSBarry Smith .keywords: SNES, nonlinear, solve 16709b94acceSBarry Smith 167170e92668SMatthew Knepley .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs(), SNESSetSolution() 16729b94acceSBarry Smith @*/ 1673f69a0ea3SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 16749b94acceSBarry Smith { 1675dfbe8321SBarry Smith PetscErrorCode ierr; 1676f1af5d2fSBarry Smith PetscTruth flg; 1677eabae89aSBarry Smith char filename[PETSC_MAX_PATH_LEN]; 1678eabae89aSBarry Smith PetscViewer viewer; 1679052efed2SBarry Smith 16803a40ed3dSBarry Smith PetscFunctionBegin; 16814482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16821302d50aSBarry Smith if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 168374637425SBarry Smith 1684f69a0ea3SMatthew Knepley if (b) { 1685f69a0ea3SMatthew Knepley ierr = SNESSetRhs(snes, b); CHKERRQ(ierr); 16861096aae1SMatthew Knepley if (!snes->vec_func) { 16871096aae1SMatthew Knepley Vec r; 16881096aae1SMatthew Knepley 16891096aae1SMatthew Knepley ierr = VecDuplicate(b, &r); CHKERRQ(ierr); 16901096aae1SMatthew Knepley ierr = SNESSetFunction(snes, r, PETSC_NULL, PETSC_NULL); CHKERRQ(ierr); 16911096aae1SMatthew Knepley } 1692f69a0ea3SMatthew Knepley } 169370e92668SMatthew Knepley if (x) { 1694f69a0ea3SMatthew Knepley PetscValidHeaderSpecific(x,VEC_COOKIE,3); 1695f69a0ea3SMatthew Knepley PetscCheckSameComm(snes,1,x,3); 169670e92668SMatthew Knepley } else { 169770e92668SMatthew Knepley ierr = SNESGetSolution(snes, &x); CHKERRQ(ierr); 169870e92668SMatthew Knepley if (!x) { 169970e92668SMatthew Knepley ierr = VecDuplicate(snes->vec_func_always, &x); CHKERRQ(ierr); 170070e92668SMatthew Knepley } 170170e92668SMatthew Knepley } 170270e92668SMatthew Knepley snes->vec_sol = snes->vec_sol_always = x; 170370e92668SMatthew Knepley if (!snes->setupcalled) { 170470e92668SMatthew Knepley ierr = SNESSetUp(snes);CHKERRQ(ierr); 170570e92668SMatthew Knepley } 1706abc0a331SBarry Smith if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1707d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 170850ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1709d5e45103SBarry Smith 1710d5e45103SBarry Smith ierr = PetscExceptionTry1((*(snes)->solve)(snes),PETSC_ERR_ARG_DOMAIN); 1711d5e45103SBarry Smith if (PetscExceptionValue(ierr)) { 171238f152feSBarry Smith /* this means that a caller above me has also tryed this exception so I don't handle it here, pass it up */ 171319717074SBarry Smith PetscErrorCode pierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(pierr); 1714d5e45103SBarry Smith } else if (PetscExceptionCaught(ierr,PETSC_ERR_ARG_DOMAIN)) { 1715d5e45103SBarry Smith /* translate exception into SNES not converged reason */ 1716d5e45103SBarry Smith snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 171738f152feSBarry Smith ierr = 0; 171838f152feSBarry Smith } 171938f152feSBarry Smith CHKERRQ(ierr); 1720d5e45103SBarry Smith 1721d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1722eabae89aSBarry Smith ierr = PetscOptionsGetString(snes->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 1723eabae89aSBarry Smith if (flg && !PetscPreLoadingOn) { 1724eabae89aSBarry Smith ierr = PetscViewerASCIIOpen(snes->comm,filename,&viewer);CHKERRQ(ierr); 1725eabae89aSBarry Smith ierr = SNESView(snes,viewer);CHKERRQ(ierr); 1726eabae89aSBarry Smith ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 1727eabae89aSBarry Smith } 1728eabae89aSBarry Smith 1729da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1730da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 17315968eb51SBarry Smith if (snes->printreason) { 17325968eb51SBarry Smith if (snes->reason > 0) { 17339dcbbd2bSBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 17345968eb51SBarry Smith } else { 17359dcbbd2bSBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 17365968eb51SBarry Smith } 17375968eb51SBarry Smith } 17385968eb51SBarry Smith 17393a40ed3dSBarry Smith PetscFunctionReturn(0); 17409b94acceSBarry Smith } 17419b94acceSBarry Smith 17429b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 17439b94acceSBarry Smith 17444a2ae208SSatish Balay #undef __FUNCT__ 17454a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 174682bf6240SBarry Smith /*@C 17474b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 17489b94acceSBarry Smith 1749fee21e36SBarry Smith Collective on SNES 1750fee21e36SBarry Smith 1751c7afd0dbSLois Curfman McInnes Input Parameters: 1752c7afd0dbSLois Curfman McInnes + snes - the SNES context 1753454a90a3SBarry Smith - type - a known method 1754c7afd0dbSLois Curfman McInnes 1755c7afd0dbSLois Curfman McInnes Options Database Key: 1756454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1757c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1758ae12b187SLois Curfman McInnes 17599b94acceSBarry Smith Notes: 1760e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 17614b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1762c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17634b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1764c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17659b94acceSBarry Smith 1766ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1767ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1768ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1769ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1770ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1771ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1772ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1773ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1774ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1775b0a32e0cSBarry Smith appropriate method. 177636851e7fSLois Curfman McInnes 177736851e7fSLois Curfman McInnes Level: intermediate 1778a703fe33SLois Curfman McInnes 1779454a90a3SBarry Smith .keywords: SNES, set, type 1780435da068SBarry Smith 1781435da068SBarry Smith .seealso: SNESType, SNESCreate() 1782435da068SBarry Smith 17839b94acceSBarry Smith @*/ 1784e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,SNESType type) 17859b94acceSBarry Smith { 1786dfbe8321SBarry Smith PetscErrorCode ierr,(*r)(SNES); 17876831982aSBarry Smith PetscTruth match; 17883a40ed3dSBarry Smith 17893a40ed3dSBarry Smith PetscFunctionBegin; 17904482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17914482741eSBarry Smith PetscValidCharPointer(type,2); 179282bf6240SBarry Smith 17936831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 17940f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 179582bf6240SBarry Smith 179682bf6240SBarry Smith if (snes->setupcalled) { 17974dc4c822SBarry Smith snes->setupcalled = PETSC_FALSE; 1798e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 179982bf6240SBarry Smith snes->data = 0; 180082bf6240SBarry Smith } 18017d1a2b2bSBarry Smith 18029b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 180382bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1804b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1805958c9bccSBarry Smith if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 180605b42c5fSBarry Smith ierr = PetscFree(snes->data);CHKERRQ(ierr); 180782bf6240SBarry Smith snes->data = 0; 18083a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1809454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 18103a40ed3dSBarry Smith PetscFunctionReturn(0); 18119b94acceSBarry Smith } 18129b94acceSBarry Smith 1813a847f771SSatish Balay 18149b94acceSBarry Smith /* --------------------------------------------------------------------- */ 18154a2ae208SSatish Balay #undef __FUNCT__ 18164a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 181752baeb72SSatish Balay /*@ 18189b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1819f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 18209b94acceSBarry Smith 1821fee21e36SBarry Smith Not Collective 1822fee21e36SBarry Smith 182336851e7fSLois Curfman McInnes Level: advanced 182436851e7fSLois Curfman McInnes 18259b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 18269b94acceSBarry Smith 18279b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 18289b94acceSBarry Smith @*/ 182963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 18309b94acceSBarry Smith { 1831dfbe8321SBarry Smith PetscErrorCode ierr; 183282bf6240SBarry Smith 18333a40ed3dSBarry Smith PetscFunctionBegin; 183482bf6240SBarry Smith if (SNESList) { 1835b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 183682bf6240SBarry Smith SNESList = 0; 18379b94acceSBarry Smith } 18384c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 18393a40ed3dSBarry Smith PetscFunctionReturn(0); 18409b94acceSBarry Smith } 18419b94acceSBarry Smith 18424a2ae208SSatish Balay #undef __FUNCT__ 18434a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 18449b94acceSBarry Smith /*@C 18459a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 18469b94acceSBarry Smith 1847c7afd0dbSLois Curfman McInnes Not Collective 1848c7afd0dbSLois Curfman McInnes 18499b94acceSBarry Smith Input Parameter: 18504b0e389bSBarry Smith . snes - nonlinear solver context 18519b94acceSBarry Smith 18529b94acceSBarry Smith Output Parameter: 18533a7fca6bSBarry Smith . type - SNES method (a character string) 18549b94acceSBarry Smith 185536851e7fSLois Curfman McInnes Level: intermediate 185636851e7fSLois Curfman McInnes 1857454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 18589b94acceSBarry Smith @*/ 185963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,SNESType *type) 18609b94acceSBarry Smith { 18613a40ed3dSBarry Smith PetscFunctionBegin; 18624482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18634482741eSBarry Smith PetscValidPointer(type,2); 1864454a90a3SBarry Smith *type = snes->type_name; 18653a40ed3dSBarry Smith PetscFunctionReturn(0); 18669b94acceSBarry Smith } 18679b94acceSBarry Smith 18684a2ae208SSatish Balay #undef __FUNCT__ 18694a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 187052baeb72SSatish Balay /*@ 18719b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 18729b94acceSBarry Smith stored. 18739b94acceSBarry Smith 1874c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1875c7afd0dbSLois Curfman McInnes 18769b94acceSBarry Smith Input Parameter: 18779b94acceSBarry Smith . snes - the SNES context 18789b94acceSBarry Smith 18799b94acceSBarry Smith Output Parameter: 18809b94acceSBarry Smith . x - the solution 18819b94acceSBarry Smith 188270e92668SMatthew Knepley Level: intermediate 188336851e7fSLois Curfman McInnes 18849b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 18859b94acceSBarry Smith 188670e92668SMatthew Knepley .seealso: SNESSetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 18879b94acceSBarry Smith @*/ 188863dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 18899b94acceSBarry Smith { 18903a40ed3dSBarry Smith PetscFunctionBegin; 18914482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18924482741eSBarry Smith PetscValidPointer(x,2); 18939b94acceSBarry Smith *x = snes->vec_sol_always; 18943a40ed3dSBarry Smith PetscFunctionReturn(0); 18959b94acceSBarry Smith } 18969b94acceSBarry Smith 18974a2ae208SSatish Balay #undef __FUNCT__ 189870e92668SMatthew Knepley #define __FUNCT__ "SNESSetSolution" 18997e4bb74cSBarry Smith /*@ 190070e92668SMatthew Knepley SNESSetSolution - Sets the vector where the approximate solution is stored. 190170e92668SMatthew Knepley 190270e92668SMatthew Knepley Not Collective, but Vec is parallel if SNES is parallel 190370e92668SMatthew Knepley 190470e92668SMatthew Knepley Input Parameters: 190570e92668SMatthew Knepley + snes - the SNES context 190670e92668SMatthew Knepley - x - the solution 190770e92668SMatthew Knepley 190870e92668SMatthew Knepley Output Parameter: 190970e92668SMatthew Knepley 191070e92668SMatthew Knepley Level: intermediate 191170e92668SMatthew Knepley 191242219521SBarry Smith Notes: this is not normally used, rather one simply calls SNESSolve() with 191342219521SBarry Smith the appropriate solution vector. 191442219521SBarry Smith 191570e92668SMatthew Knepley .keywords: SNES, nonlinear, set, solution 191670e92668SMatthew Knepley 191770e92668SMatthew Knepley .seealso: SNESGetSolution(), SNESGetFunction(), SNESGetSolutionUpdate() 191870e92668SMatthew Knepley @*/ 191970e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESSetSolution(SNES snes,Vec x) 192070e92668SMatthew Knepley { 192170e92668SMatthew Knepley PetscFunctionBegin; 192270e92668SMatthew Knepley PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 192370e92668SMatthew Knepley PetscValidHeaderSpecific(x,VEC_COOKIE,2); 192470e92668SMatthew Knepley PetscCheckSameComm(snes,1,x,2); 192570e92668SMatthew Knepley snes->vec_sol_always = x; 192670e92668SMatthew Knepley PetscFunctionReturn(0); 192770e92668SMatthew Knepley } 192870e92668SMatthew Knepley 192970e92668SMatthew Knepley #undef __FUNCT__ 19304a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 193152baeb72SSatish Balay /*@ 19329b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 19339b94acceSBarry Smith stored. 19349b94acceSBarry Smith 1935c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1936c7afd0dbSLois Curfman McInnes 19379b94acceSBarry Smith Input Parameter: 19389b94acceSBarry Smith . snes - the SNES context 19399b94acceSBarry Smith 19409b94acceSBarry Smith Output Parameter: 19419b94acceSBarry Smith . x - the solution update 19429b94acceSBarry Smith 194336851e7fSLois Curfman McInnes Level: advanced 194436851e7fSLois Curfman McInnes 19459b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 19469b94acceSBarry Smith 19479b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 19489b94acceSBarry Smith @*/ 194963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 19509b94acceSBarry Smith { 19513a40ed3dSBarry Smith PetscFunctionBegin; 19524482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19534482741eSBarry Smith PetscValidPointer(x,2); 19549b94acceSBarry Smith *x = snes->vec_sol_update_always; 19553a40ed3dSBarry Smith PetscFunctionReturn(0); 19569b94acceSBarry Smith } 19579b94acceSBarry Smith 19584a2ae208SSatish Balay #undef __FUNCT__ 19594a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 19609b94acceSBarry Smith /*@C 19613638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 19629b94acceSBarry Smith 1963c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1964c7afd0dbSLois Curfman McInnes 19659b94acceSBarry Smith Input Parameter: 19669b94acceSBarry Smith . snes - the SNES context 19679b94acceSBarry Smith 19689b94acceSBarry Smith Output Parameter: 19697bf4e008SBarry Smith + r - the function (or PETSC_NULL) 197070e92668SMatthew Knepley . func - the function (or PETSC_NULL) 197170e92668SMatthew Knepley - ctx - the function context (or PETSC_NULL) 19729b94acceSBarry Smith 197336851e7fSLois Curfman McInnes Level: advanced 197436851e7fSLois Curfman McInnes 1975a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 19769b94acceSBarry Smith 19774b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 19789b94acceSBarry Smith @*/ 197970e92668SMatthew Knepley PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 19809b94acceSBarry Smith { 19813a40ed3dSBarry Smith PetscFunctionBegin; 19824482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19837bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 198400036973SBarry Smith if (func) *func = snes->computefunction; 198570e92668SMatthew Knepley if (ctx) *ctx = snes->funP; 19863a40ed3dSBarry Smith PetscFunctionReturn(0); 19879b94acceSBarry Smith } 19889b94acceSBarry Smith 19894a2ae208SSatish Balay #undef __FUNCT__ 19904a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 19913c7409f5SSatish Balay /*@C 19923c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1993d850072dSLois Curfman McInnes SNES options in the database. 19943c7409f5SSatish Balay 1995fee21e36SBarry Smith Collective on SNES 1996fee21e36SBarry Smith 1997c7afd0dbSLois Curfman McInnes Input Parameter: 1998c7afd0dbSLois Curfman McInnes + snes - the SNES context 1999c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2000c7afd0dbSLois Curfman McInnes 2001d850072dSLois Curfman McInnes Notes: 2002a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2003c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2004d850072dSLois Curfman McInnes 200536851e7fSLois Curfman McInnes Level: advanced 200636851e7fSLois Curfman McInnes 20073c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 2008a86d99e1SLois Curfman McInnes 2009a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 20103c7409f5SSatish Balay @*/ 201163dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 20123c7409f5SSatish Balay { 2013dfbe8321SBarry Smith PetscErrorCode ierr; 20143c7409f5SSatish Balay 20153a40ed3dSBarry Smith PetscFunctionBegin; 20164482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2017639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 201894b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20193a40ed3dSBarry Smith PetscFunctionReturn(0); 20203c7409f5SSatish Balay } 20213c7409f5SSatish Balay 20224a2ae208SSatish Balay #undef __FUNCT__ 20234a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 20243c7409f5SSatish Balay /*@C 2025f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2026d850072dSLois Curfman McInnes SNES options in the database. 20273c7409f5SSatish Balay 2028fee21e36SBarry Smith Collective on SNES 2029fee21e36SBarry Smith 2030c7afd0dbSLois Curfman McInnes Input Parameters: 2031c7afd0dbSLois Curfman McInnes + snes - the SNES context 2032c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2033c7afd0dbSLois Curfman McInnes 2034d850072dSLois Curfman McInnes Notes: 2035a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2036c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2037d850072dSLois Curfman McInnes 203836851e7fSLois Curfman McInnes Level: advanced 203936851e7fSLois Curfman McInnes 20403c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2041a86d99e1SLois Curfman McInnes 2042a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20433c7409f5SSatish Balay @*/ 204463dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 20453c7409f5SSatish Balay { 2046dfbe8321SBarry Smith PetscErrorCode ierr; 20473c7409f5SSatish Balay 20483a40ed3dSBarry Smith PetscFunctionBegin; 20494482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2050639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 205194b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20523a40ed3dSBarry Smith PetscFunctionReturn(0); 20533c7409f5SSatish Balay } 20543c7409f5SSatish Balay 20554a2ae208SSatish Balay #undef __FUNCT__ 20564a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 20579ab63eb5SSatish Balay /*@C 20583c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20593c7409f5SSatish Balay SNES options in the database. 20603c7409f5SSatish Balay 2061c7afd0dbSLois Curfman McInnes Not Collective 2062c7afd0dbSLois Curfman McInnes 20633c7409f5SSatish Balay Input Parameter: 20643c7409f5SSatish Balay . snes - the SNES context 20653c7409f5SSatish Balay 20663c7409f5SSatish Balay Output Parameter: 20673c7409f5SSatish Balay . prefix - pointer to the prefix string used 20683c7409f5SSatish Balay 20699ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 20709ab63eb5SSatish Balay sufficient length to hold the prefix. 20719ab63eb5SSatish Balay 207236851e7fSLois Curfman McInnes Level: advanced 207336851e7fSLois Curfman McInnes 20743c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2075a86d99e1SLois Curfman McInnes 2076a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20773c7409f5SSatish Balay @*/ 2078e060cb09SBarry Smith PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 20793c7409f5SSatish Balay { 2080dfbe8321SBarry Smith PetscErrorCode ierr; 20813c7409f5SSatish Balay 20823a40ed3dSBarry Smith PetscFunctionBegin; 20834482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2084639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 20853a40ed3dSBarry Smith PetscFunctionReturn(0); 20863c7409f5SSatish Balay } 20873c7409f5SSatish Balay 2088b2002411SLois Curfman McInnes 20894a2ae208SSatish Balay #undef __FUNCT__ 20904a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 20913cea93caSBarry Smith /*@C 20923cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 20933cea93caSBarry Smith 20947f6c08e0SMatthew Knepley Level: advanced 20953cea93caSBarry Smith @*/ 209663dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2097b2002411SLois Curfman McInnes { 2098e2d1d2b7SBarry Smith char fullname[PETSC_MAX_PATH_LEN]; 2099dfbe8321SBarry Smith PetscErrorCode ierr; 2100b2002411SLois Curfman McInnes 2101b2002411SLois Curfman McInnes PetscFunctionBegin; 2102b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2103c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2104b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2105b2002411SLois Curfman McInnes } 2106da9b6338SBarry Smith 2107da9b6338SBarry Smith #undef __FUNCT__ 2108da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 210963dd3a1aSKris Buschelman PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2110da9b6338SBarry Smith { 2111dfbe8321SBarry Smith PetscErrorCode ierr; 211277431f27SBarry Smith PetscInt N,i,j; 2113da9b6338SBarry Smith Vec u,uh,fh; 2114da9b6338SBarry Smith PetscScalar value; 2115da9b6338SBarry Smith PetscReal norm; 2116da9b6338SBarry Smith 2117da9b6338SBarry Smith PetscFunctionBegin; 2118da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2119da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2120da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2121da9b6338SBarry Smith 2122da9b6338SBarry Smith /* currently only works for sequential */ 2123da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2124da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2125da9b6338SBarry Smith for (i=0; i<N; i++) { 2126da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 212777431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2128da9b6338SBarry Smith for (j=-10; j<11; j++) { 2129ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2130da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 21313ab0aad5SBarry Smith ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2132da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 213377431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2134da9b6338SBarry Smith value = -value; 2135da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2136da9b6338SBarry Smith } 2137da9b6338SBarry Smith } 2138da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2139da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2140da9b6338SBarry Smith PetscFunctionReturn(0); 2141da9b6338SBarry Smith } 2142