19b94acceSBarry Smith 2e090d566SSatish Balay #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 39b94acceSBarry Smith 44c49b128SBarry Smith PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 58ba1e511SMatthew Knepley PetscFList SNESList = PETSC_NULL; 68ba1e511SMatthew Knepley 78ba1e511SMatthew Knepley /* Logging support */ 86849ba73SBarry Smith PetscCookie SNES_COOKIE = 0; 96849ba73SBarry Smith PetscEvent SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0; 10a09944afSBarry Smith 11a09944afSBarry Smith #undef __FUNCT__ 124a2ae208SSatish Balay #define __FUNCT__ "SNESView" 137e2c5f70SBarry Smith /*@C 149b94acceSBarry Smith SNESView - Prints the SNES data structure. 159b94acceSBarry Smith 164c49b128SBarry Smith Collective on SNES 17fee21e36SBarry Smith 18c7afd0dbSLois Curfman McInnes Input Parameters: 19c7afd0dbSLois Curfman McInnes + SNES - the SNES context 20c7afd0dbSLois Curfman McInnes - viewer - visualization context 21c7afd0dbSLois Curfman McInnes 229b94acceSBarry Smith Options Database Key: 23c8a8ba5cSLois Curfman McInnes . -snes_view - Calls SNESView() at end of SNESSolve() 249b94acceSBarry Smith 259b94acceSBarry Smith Notes: 269b94acceSBarry Smith The available visualization contexts include 27b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 28b0a32e0cSBarry Smith - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 29c8a8ba5cSLois Curfman McInnes output where only the first processor opens 30c8a8ba5cSLois Curfman McInnes the file. All other processors send their 31c8a8ba5cSLois Curfman McInnes data to the first processor to print. 329b94acceSBarry Smith 333e081fefSLois Curfman McInnes The user can open an alternative visualization context with 34b0a32e0cSBarry Smith PetscViewerASCIIOpen() - output to a specified file. 359b94acceSBarry Smith 3636851e7fSLois Curfman McInnes Level: beginner 3736851e7fSLois Curfman McInnes 389b94acceSBarry Smith .keywords: SNES, view 399b94acceSBarry Smith 40b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen() 419b94acceSBarry Smith @*/ 42dfbe8321SBarry Smith PetscErrorCode SNESView(SNES snes,PetscViewer viewer) 439b94acceSBarry Smith { 449b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 45dfbe8321SBarry Smith PetscErrorCode ierr; 4694b7f48cSBarry Smith KSP ksp; 47454a90a3SBarry Smith char *type; 4832077d6dSBarry Smith PetscTruth iascii,isstring; 499b94acceSBarry Smith 503a40ed3dSBarry Smith PetscFunctionBegin; 514482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 52b0a32e0cSBarry Smith if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm); 534482741eSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 54c9780b6fSBarry Smith PetscCheckSameComm(snes,1,viewer,2); 5574679c65SBarry Smith 5632077d6dSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 57b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 5832077d6dSBarry Smith if (iascii) { 593a7fca6bSBarry Smith if (snes->prefix) { 603a7fca6bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr); 613a7fca6bSBarry Smith } else { 62b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 633a7fca6bSBarry Smith } 64454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 65454a90a3SBarry Smith if (type) { 66b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 67184914b5SBarry Smith } else { 68b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 69184914b5SBarry Smith } 700ef38995SBarry Smith if (snes->view) { 71b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 720ef38995SBarry Smith ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 73b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 740ef38995SBarry Smith } 7577431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 7677d8c4bbSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n", 7770441072SBarry Smith snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr); 7877431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 7977431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 809b94acceSBarry Smith if (snes->ksp_ewconv) { 819b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 829b94acceSBarry Smith if (kctx) { 8377431f27SBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 84b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 85b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 869b94acceSBarry Smith } 879b94acceSBarry Smith } 880f5bd95cSBarry Smith } else if (isstring) { 89454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 90b0a32e0cSBarry Smith ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 9119bcc07fSBarry Smith } 9294b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 93b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 9494b7f48cSBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 95b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 963a40ed3dSBarry Smith PetscFunctionReturn(0); 979b94acceSBarry Smith } 989b94acceSBarry Smith 9976b2cf59SMatthew Knepley /* 10076b2cf59SMatthew Knepley We retain a list of functions that also take SNES command 10176b2cf59SMatthew Knepley line options. These are called at the end SNESSetFromOptions() 10276b2cf59SMatthew Knepley */ 10376b2cf59SMatthew Knepley #define MAXSETFROMOPTIONS 5 104a7cc72afSBarry Smith static PetscInt numberofsetfromoptions; 1056849ba73SBarry Smith static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 10676b2cf59SMatthew Knepley 107e74ef692SMatthew Knepley #undef __FUNCT__ 108e74ef692SMatthew Knepley #define __FUNCT__ "SNESAddOptionsChecker" 109ac226902SBarry Smith /*@C 11076b2cf59SMatthew Knepley SNESAddOptionsChecker - Adds an additional function to check for SNES options. 11176b2cf59SMatthew Knepley 11276b2cf59SMatthew Knepley Not Collective 11376b2cf59SMatthew Knepley 11476b2cf59SMatthew Knepley Input Parameter: 11576b2cf59SMatthew Knepley . snescheck - function that checks for options 11676b2cf59SMatthew Knepley 11776b2cf59SMatthew Knepley Level: developer 11876b2cf59SMatthew Knepley 11976b2cf59SMatthew Knepley .seealso: SNESSetFromOptions() 12076b2cf59SMatthew Knepley @*/ 1216849ba73SBarry Smith PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 12276b2cf59SMatthew Knepley { 12376b2cf59SMatthew Knepley PetscFunctionBegin; 12476b2cf59SMatthew Knepley if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 12577431f27SBarry Smith SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 12676b2cf59SMatthew Knepley } 12776b2cf59SMatthew Knepley othersetfromoptions[numberofsetfromoptions++] = snescheck; 12876b2cf59SMatthew Knepley PetscFunctionReturn(0); 12976b2cf59SMatthew Knepley } 13076b2cf59SMatthew Knepley 1314a2ae208SSatish Balay #undef __FUNCT__ 1324a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions" 1339b94acceSBarry Smith /*@ 13494b7f48cSBarry Smith SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 1359b94acceSBarry Smith 136c7afd0dbSLois Curfman McInnes Collective on SNES 137c7afd0dbSLois Curfman McInnes 1389b94acceSBarry Smith Input Parameter: 1399b94acceSBarry Smith . snes - the SNES context 1409b94acceSBarry Smith 14136851e7fSLois Curfman McInnes Options Database Keys: 1426831982aSBarry Smith + -snes_type <type> - ls, tr, umls, umtr, test 14382738288SBarry Smith . -snes_stol - convergence tolerance in terms of the norm 14482738288SBarry Smith of the change in the solution between steps 14570441072SBarry Smith . -snes_atol <abstol> - absolute tolerance of residual norm 146b39c3a46SLois Curfman McInnes . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 147b39c3a46SLois Curfman McInnes . -snes_max_it <max_it> - maximum number of iterations 148b39c3a46SLois Curfman McInnes . -snes_max_funcs <max_funcs> - maximum number of function evaluations 14950ffb88aSMatthew Knepley . -snes_max_fail <max_fail> - maximum number of failures 150b39c3a46SLois Curfman McInnes . -snes_trtol <trtol> - trust region tolerance 1512492ecdbSBarry Smith . -snes_no_convergence_test - skip convergence test in nonlinear 15282738288SBarry Smith solver; hence iterations will continue until max_it 1531fbbfb26SLois Curfman McInnes or some other criterion is reached. Saves expense 15482738288SBarry Smith of convergence test 15582738288SBarry Smith . -snes_monitor - prints residual norm at each iteration 156d132466eSBarry Smith . -snes_vecmonitor - plots solution at each iteration 157d132466eSBarry Smith . -snes_vecmonitor_update - plots update to solution at each iteration 15882738288SBarry Smith . -snes_xmonitor - plots residual norm at each iteration 159e24b481bSBarry Smith . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 1605968eb51SBarry Smith . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 1615968eb51SBarry Smith - -snes_print_converged_reason - print the reason for convergence/divergence after each solve 16282738288SBarry Smith 16382738288SBarry Smith Options Database for Eisenstat-Walker method: 1644b27c08aSLois Curfman McInnes + -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence 1654b27c08aSLois Curfman McInnes . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 16636851e7fSLois Curfman McInnes . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 16736851e7fSLois Curfman McInnes . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 16836851e7fSLois Curfman McInnes . -snes_ksp_ew_gamma <gamma> - Sets gamma 16936851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha <alpha> - Sets alpha 17036851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 17136851e7fSLois Curfman McInnes - -snes_ksp_ew_threshold <threshold> - Sets threshold 17282738288SBarry Smith 17311ca99fdSLois Curfman McInnes Notes: 17411ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 17511ca99fdSLois Curfman McInnes the users manual. 17683e2fdc7SBarry Smith 17736851e7fSLois Curfman McInnes Level: beginner 17836851e7fSLois Curfman McInnes 1799b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1809b94acceSBarry Smith 18169ed319cSSatish Balay .seealso: SNESSetOptionsPrefix() 1829b94acceSBarry Smith @*/ 183dfbe8321SBarry Smith PetscErrorCode SNESSetFromOptions(SNES snes) 1849b94acceSBarry Smith { 18594b7f48cSBarry Smith KSP ksp; 186186905e3SBarry Smith SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 187f1af5d2fSBarry Smith PetscTruth flg; 188dfbe8321SBarry Smith PetscErrorCode ierr; 18977431f27SBarry Smith PetscInt i; 1902fc52814SBarry Smith const char *deft; 1912fc52814SBarry Smith char type[256]; 1929b94acceSBarry Smith 1933a40ed3dSBarry Smith PetscFunctionBegin; 1944482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 195ca161407SBarry Smith 196b0a32e0cSBarry Smith ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 197186905e3SBarry Smith if (snes->type_name) { 198186905e3SBarry Smith deft = snes->type_name; 199186905e3SBarry Smith } else { 2004b27c08aSLois Curfman McInnes deft = SNESLS; 201d64ed03dSBarry Smith } 2024bbc92c1SBarry Smith 203186905e3SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 204b0a32e0cSBarry Smith ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 205d64ed03dSBarry Smith if (flg) { 206186905e3SBarry Smith ierr = SNESSetType(snes,type);CHKERRQ(ierr); 207186905e3SBarry Smith } else if (!snes->type_name) { 208186905e3SBarry Smith ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 209d64ed03dSBarry Smith } 210909c8a9fSBarry Smith ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 21193c39befSBarry Smith 21287828ca2SBarry Smith ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 21370441072SBarry Smith ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr); 214186905e3SBarry Smith 21587828ca2SBarry Smith ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 216b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 217b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 21850ffb88aSMatthew Knepley ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 2195968eb51SBarry Smith ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr); 2205968eb51SBarry Smith if (flg) { 2215968eb51SBarry Smith snes->printreason = PETSC_TRUE; 2225968eb51SBarry Smith } 223186905e3SBarry Smith 224b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr); 225186905e3SBarry Smith 226b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 22787828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 22887828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 22987828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 23087828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 23187828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 23287828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 233186905e3SBarry Smith 234b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr); 23593c39befSBarry Smith if (flg) {snes->converged = 0;} 236b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr); 2375cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 238b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr); 239b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);} 2403a7fca6bSBarry Smith ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr); 2413a7fca6bSBarry Smith if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);} 242b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr); 243b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);} 244b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr); 245b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 246b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr); 2477c922b88SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);} 2485ed2d596SBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr); 2495ed2d596SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);} 250b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr); 251186905e3SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 252e24b481bSBarry Smith 253b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 2544b27c08aSLois Curfman McInnes if (flg) { 255186905e3SBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 256b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 2579b94acceSBarry Smith } 258639f9d9dSBarry Smith 25976b2cf59SMatthew Knepley for(i = 0; i < numberofsetfromoptions; i++) { 26076b2cf59SMatthew Knepley ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 26176b2cf59SMatthew Knepley } 26276b2cf59SMatthew Knepley 263186905e3SBarry Smith if (snes->setfromoptions) { 264186905e3SBarry Smith ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 265639f9d9dSBarry Smith } 266639f9d9dSBarry Smith 267b0a32e0cSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 2684bbc92c1SBarry Smith 26994b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 27094b7f48cSBarry Smith ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 27193993e2dSLois Curfman McInnes 2723a40ed3dSBarry Smith PetscFunctionReturn(0); 2739b94acceSBarry Smith } 2749b94acceSBarry Smith 275a847f771SSatish Balay 2764a2ae208SSatish Balay #undef __FUNCT__ 2774a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext" 2789b94acceSBarry Smith /*@ 2799b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 2809b94acceSBarry Smith the nonlinear solvers. 2819b94acceSBarry Smith 282fee21e36SBarry Smith Collective on SNES 283fee21e36SBarry Smith 284c7afd0dbSLois Curfman McInnes Input Parameters: 285c7afd0dbSLois Curfman McInnes + snes - the SNES context 286c7afd0dbSLois Curfman McInnes - usrP - optional user context 287c7afd0dbSLois Curfman McInnes 28836851e7fSLois Curfman McInnes Level: intermediate 28936851e7fSLois Curfman McInnes 2909b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 2919b94acceSBarry Smith 2929b94acceSBarry Smith .seealso: SNESGetApplicationContext() 2939b94acceSBarry Smith @*/ 294dfbe8321SBarry Smith PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP) 2959b94acceSBarry Smith { 2963a40ed3dSBarry Smith PetscFunctionBegin; 2974482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2989b94acceSBarry Smith snes->user = usrP; 2993a40ed3dSBarry Smith PetscFunctionReturn(0); 3009b94acceSBarry Smith } 30174679c65SBarry Smith 3024a2ae208SSatish Balay #undef __FUNCT__ 3034a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext" 3049b94acceSBarry Smith /*@C 3059b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3069b94acceSBarry Smith nonlinear solvers. 3079b94acceSBarry Smith 308c7afd0dbSLois Curfman McInnes Not Collective 309c7afd0dbSLois Curfman McInnes 3109b94acceSBarry Smith Input Parameter: 3119b94acceSBarry Smith . snes - SNES context 3129b94acceSBarry Smith 3139b94acceSBarry Smith Output Parameter: 3149b94acceSBarry Smith . usrP - user context 3159b94acceSBarry Smith 31636851e7fSLois Curfman McInnes Level: intermediate 31736851e7fSLois Curfman McInnes 3189b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3199b94acceSBarry Smith 3209b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3219b94acceSBarry Smith @*/ 322dfbe8321SBarry Smith PetscErrorCode SNESGetApplicationContext(SNES snes,void **usrP) 3239b94acceSBarry Smith { 3243a40ed3dSBarry Smith PetscFunctionBegin; 3254482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 3269b94acceSBarry Smith *usrP = snes->user; 3273a40ed3dSBarry Smith PetscFunctionReturn(0); 3289b94acceSBarry Smith } 32974679c65SBarry Smith 3304a2ae208SSatish Balay #undef __FUNCT__ 3314a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber" 3329b94acceSBarry Smith /*@ 333c8228a4eSBarry Smith SNESGetIterationNumber - Gets the number of nonlinear iterations completed 334c8228a4eSBarry Smith at this time. 3359b94acceSBarry Smith 336c7afd0dbSLois Curfman McInnes Not Collective 337c7afd0dbSLois Curfman McInnes 3389b94acceSBarry Smith Input Parameter: 3399b94acceSBarry Smith . snes - SNES context 3409b94acceSBarry Smith 3419b94acceSBarry Smith Output Parameter: 3429b94acceSBarry Smith . iter - iteration number 3439b94acceSBarry Smith 344c8228a4eSBarry Smith Notes: 345c8228a4eSBarry Smith For example, during the computation of iteration 2 this would return 1. 346c8228a4eSBarry Smith 347c8228a4eSBarry Smith This is useful for using lagged Jacobians (where one does not recompute the 34808405cd6SLois Curfman McInnes Jacobian at each SNES iteration). For example, the code 34908405cd6SLois Curfman McInnes .vb 35008405cd6SLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&it); 35108405cd6SLois Curfman McInnes if (!(it % 2)) { 35208405cd6SLois Curfman McInnes [compute Jacobian here] 35308405cd6SLois Curfman McInnes } 35408405cd6SLois Curfman McInnes .ve 355c8228a4eSBarry Smith can be used in your ComputeJacobian() function to cause the Jacobian to be 35608405cd6SLois Curfman McInnes recomputed every second SNES iteration. 357c8228a4eSBarry Smith 35836851e7fSLois Curfman McInnes Level: intermediate 35936851e7fSLois Curfman McInnes 3609b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3619b94acceSBarry Smith @*/ 36277431f27SBarry Smith PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt* iter) 3639b94acceSBarry Smith { 3643a40ed3dSBarry Smith PetscFunctionBegin; 3654482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 3664482741eSBarry Smith PetscValidIntPointer(iter,2); 3679b94acceSBarry Smith *iter = snes->iter; 3683a40ed3dSBarry Smith PetscFunctionReturn(0); 3699b94acceSBarry Smith } 37074679c65SBarry Smith 3714a2ae208SSatish Balay #undef __FUNCT__ 3724a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm" 3739b94acceSBarry Smith /*@ 3749b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3759b94acceSBarry Smith with SNESSSetFunction(). 3769b94acceSBarry Smith 377c7afd0dbSLois Curfman McInnes Collective on SNES 378c7afd0dbSLois Curfman McInnes 3799b94acceSBarry Smith Input Parameter: 3809b94acceSBarry Smith . snes - SNES context 3819b94acceSBarry Smith 3829b94acceSBarry Smith Output Parameter: 3839b94acceSBarry Smith . fnorm - 2-norm of function 3849b94acceSBarry Smith 38536851e7fSLois Curfman McInnes Level: intermediate 38636851e7fSLois Curfman McInnes 3879b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 388a86d99e1SLois Curfman McInnes 38908405cd6SLois Curfman McInnes .seealso: SNESGetFunction() 3909b94acceSBarry Smith @*/ 391dfbe8321SBarry Smith PetscErrorCode SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm) 3929b94acceSBarry Smith { 3933a40ed3dSBarry Smith PetscFunctionBegin; 3944482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 3954482741eSBarry Smith PetscValidScalarPointer(fnorm,2); 3969b94acceSBarry Smith *fnorm = snes->norm; 3973a40ed3dSBarry Smith PetscFunctionReturn(0); 3989b94acceSBarry Smith } 39974679c65SBarry Smith 4004a2ae208SSatish Balay #undef __FUNCT__ 4014a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps" 4029b94acceSBarry Smith /*@ 4039b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4049b94acceSBarry Smith attempted by the nonlinear solver. 4059b94acceSBarry Smith 406c7afd0dbSLois Curfman McInnes Not Collective 407c7afd0dbSLois Curfman McInnes 4089b94acceSBarry Smith Input Parameter: 4099b94acceSBarry Smith . snes - SNES context 4109b94acceSBarry Smith 4119b94acceSBarry Smith Output Parameter: 4129b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4139b94acceSBarry Smith 414c96a6f78SLois Curfman McInnes Notes: 415c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 416c96a6f78SLois Curfman McInnes 41736851e7fSLois Curfman McInnes Level: intermediate 41836851e7fSLois Curfman McInnes 4199b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4209b94acceSBarry Smith @*/ 42177431f27SBarry Smith PetscErrorCode SNESGetNumberUnsuccessfulSteps(SNES snes,PetscInt* nfails) 4229b94acceSBarry Smith { 4233a40ed3dSBarry Smith PetscFunctionBegin; 4244482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 4254482741eSBarry Smith PetscValidIntPointer(nfails,2); 42650ffb88aSMatthew Knepley *nfails = snes->numFailures; 42750ffb88aSMatthew Knepley PetscFunctionReturn(0); 42850ffb88aSMatthew Knepley } 42950ffb88aSMatthew Knepley 43050ffb88aSMatthew Knepley #undef __FUNCT__ 43150ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps" 43250ffb88aSMatthew Knepley /*@ 43350ffb88aSMatthew Knepley SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps 43450ffb88aSMatthew Knepley attempted by the nonlinear solver before it gives up. 43550ffb88aSMatthew Knepley 43650ffb88aSMatthew Knepley Not Collective 43750ffb88aSMatthew Knepley 43850ffb88aSMatthew Knepley Input Parameters: 43950ffb88aSMatthew Knepley + snes - SNES context 44050ffb88aSMatthew Knepley - maxFails - maximum of unsuccessful steps 44150ffb88aSMatthew Knepley 44250ffb88aSMatthew Knepley Level: intermediate 44350ffb88aSMatthew Knepley 44450ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 44550ffb88aSMatthew Knepley @*/ 44677431f27SBarry Smith PetscErrorCode SNESSetMaximumUnsuccessfulSteps(SNES snes, PetscInt maxFails) 44750ffb88aSMatthew Knepley { 44850ffb88aSMatthew Knepley PetscFunctionBegin; 4494482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 45050ffb88aSMatthew Knepley snes->maxFailures = maxFails; 45150ffb88aSMatthew Knepley PetscFunctionReturn(0); 45250ffb88aSMatthew Knepley } 45350ffb88aSMatthew Knepley 45450ffb88aSMatthew Knepley #undef __FUNCT__ 45550ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps" 45650ffb88aSMatthew Knepley /*@ 45750ffb88aSMatthew Knepley SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps 45850ffb88aSMatthew Knepley attempted by the nonlinear solver before it gives up. 45950ffb88aSMatthew Knepley 46050ffb88aSMatthew Knepley Not Collective 46150ffb88aSMatthew Knepley 46250ffb88aSMatthew Knepley Input Parameter: 46350ffb88aSMatthew Knepley . snes - SNES context 46450ffb88aSMatthew Knepley 46550ffb88aSMatthew Knepley Output Parameter: 46650ffb88aSMatthew Knepley . maxFails - maximum of unsuccessful steps 46750ffb88aSMatthew Knepley 46850ffb88aSMatthew Knepley Level: intermediate 46950ffb88aSMatthew Knepley 47050ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 47150ffb88aSMatthew Knepley @*/ 472a7cc72afSBarry Smith PetscErrorCode SNESGetMaximumUnsuccessfulSteps(SNES snes, PetscInt *maxFails) 47350ffb88aSMatthew Knepley { 47450ffb88aSMatthew Knepley PetscFunctionBegin; 4754482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 4764482741eSBarry Smith PetscValidIntPointer(maxFails,2); 47750ffb88aSMatthew Knepley *maxFails = snes->maxFailures; 4783a40ed3dSBarry Smith PetscFunctionReturn(0); 4799b94acceSBarry Smith } 480a847f771SSatish Balay 4814a2ae208SSatish Balay #undef __FUNCT__ 4824a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations" 483c96a6f78SLois Curfman McInnes /*@ 484c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 485c96a6f78SLois Curfman McInnes used by the nonlinear solver. 486c96a6f78SLois Curfman McInnes 487c7afd0dbSLois Curfman McInnes Not Collective 488c7afd0dbSLois Curfman McInnes 489c96a6f78SLois Curfman McInnes Input Parameter: 490c96a6f78SLois Curfman McInnes . snes - SNES context 491c96a6f78SLois Curfman McInnes 492c96a6f78SLois Curfman McInnes Output Parameter: 493c96a6f78SLois Curfman McInnes . lits - number of linear iterations 494c96a6f78SLois Curfman McInnes 495c96a6f78SLois Curfman McInnes Notes: 496c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 497c96a6f78SLois Curfman McInnes 49836851e7fSLois Curfman McInnes Level: intermediate 49936851e7fSLois Curfman McInnes 500c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 501c96a6f78SLois Curfman McInnes @*/ 502a7cc72afSBarry Smith PetscErrorCode SNESGetNumberLinearIterations(SNES snes,PetscInt* lits) 503c96a6f78SLois Curfman McInnes { 5043a40ed3dSBarry Smith PetscFunctionBegin; 5054482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 5064482741eSBarry Smith PetscValidIntPointer(lits,2); 507c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 5083a40ed3dSBarry Smith PetscFunctionReturn(0); 509c96a6f78SLois Curfman McInnes } 510c96a6f78SLois Curfman McInnes 5114a2ae208SSatish Balay #undef __FUNCT__ 51294b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP" 5139b94acceSBarry Smith /*@C 51494b7f48cSBarry Smith SNESGetKSP - Returns the KSP context for a SNES solver. 5159b94acceSBarry Smith 51694b7f48cSBarry Smith Not Collective, but if SNES object is parallel, then KSP object is parallel 517c7afd0dbSLois Curfman McInnes 5189b94acceSBarry Smith Input Parameter: 5199b94acceSBarry Smith . snes - the SNES context 5209b94acceSBarry Smith 5219b94acceSBarry Smith Output Parameter: 52294b7f48cSBarry Smith . ksp - the KSP context 5239b94acceSBarry Smith 5249b94acceSBarry Smith Notes: 52594b7f48cSBarry Smith The user can then directly manipulate the KSP context to set various 5269b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5279b94acceSBarry Smith KSP and PC contexts as well. 5289b94acceSBarry Smith 52936851e7fSLois Curfman McInnes Level: beginner 53036851e7fSLois Curfman McInnes 53194b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context 5329b94acceSBarry Smith 53394b7f48cSBarry Smith .seealso: KSPGetPC() 5349b94acceSBarry Smith @*/ 535dfbe8321SBarry Smith PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp) 5369b94acceSBarry Smith { 5373a40ed3dSBarry Smith PetscFunctionBegin; 5384482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 5394482741eSBarry Smith PetscValidPointer(ksp,2); 54094b7f48cSBarry Smith *ksp = snes->ksp; 5413a40ed3dSBarry Smith PetscFunctionReturn(0); 5429b94acceSBarry Smith } 54382bf6240SBarry Smith 5444a2ae208SSatish Balay #undef __FUNCT__ 5454a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc" 5466849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 547e24b481bSBarry Smith { 548e24b481bSBarry Smith PetscFunctionBegin; 549e24b481bSBarry Smith PetscFunctionReturn(0); 550e24b481bSBarry Smith } 551e24b481bSBarry Smith 5529b94acceSBarry Smith /* -----------------------------------------------------------*/ 5534a2ae208SSatish Balay #undef __FUNCT__ 5544a2ae208SSatish Balay #define __FUNCT__ "SNESCreate" 5559b94acceSBarry Smith /*@C 5569b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5579b94acceSBarry Smith 558c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 559c7afd0dbSLois Curfman McInnes 560c7afd0dbSLois Curfman McInnes Input Parameters: 561c7afd0dbSLois Curfman McInnes + comm - MPI communicator 5629b94acceSBarry Smith 5639b94acceSBarry Smith Output Parameter: 5649b94acceSBarry Smith . outsnes - the new SNES context 5659b94acceSBarry Smith 566c7afd0dbSLois Curfman McInnes Options Database Keys: 567c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 568c7afd0dbSLois Curfman McInnes and no preconditioning matrix 569c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 570c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 571c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 572c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 573c1f60f51SBarry Smith 57436851e7fSLois Curfman McInnes Level: beginner 57536851e7fSLois Curfman McInnes 5769b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5779b94acceSBarry Smith 5784b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES 5799b94acceSBarry Smith @*/ 580dfbe8321SBarry Smith PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes) 5819b94acceSBarry Smith { 582dfbe8321SBarry Smith PetscErrorCode ierr; 5839b94acceSBarry Smith SNES snes; 5849b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 58537fcc0dbSBarry Smith 5863a40ed3dSBarry Smith PetscFunctionBegin; 5874482741eSBarry Smith PetscValidPointer(outsnes,1); 5888ba1e511SMatthew Knepley *outsnes = PETSC_NULL; 5898ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 5908ba1e511SMatthew Knepley ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 5918ba1e511SMatthew Knepley #endif 5928ba1e511SMatthew Knepley 59352e6d16bSBarry Smith ierr = PetscHeaderCreate(snes,_p_SNES,PetscInt,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 594e24b481bSBarry Smith snes->bops->publish = SNESPublish_Petsc; 5959b94acceSBarry Smith snes->max_its = 50; 5969750a799SBarry Smith snes->max_funcs = 10000; 5979b94acceSBarry Smith snes->norm = 0.0; 598b4874afaSBarry Smith snes->rtol = 1.e-8; 599b4874afaSBarry Smith snes->ttol = 0.0; 60070441072SBarry Smith snes->abstol = 1.e-50; 6019b94acceSBarry Smith snes->xtol = 1.e-8; 6024b27c08aSLois Curfman McInnes snes->deltatol = 1.e-12; 6039b94acceSBarry Smith snes->nfuncs = 0; 60450ffb88aSMatthew Knepley snes->numFailures = 0; 60550ffb88aSMatthew Knepley snes->maxFailures = 1; 6067a00f4a9SLois Curfman McInnes snes->linear_its = 0; 607639f9d9dSBarry Smith snes->numbermonitors = 0; 6089b94acceSBarry Smith snes->data = 0; 6099b94acceSBarry Smith snes->view = 0; 61082bf6240SBarry Smith snes->setupcalled = 0; 611186905e3SBarry Smith snes->ksp_ewconv = PETSC_FALSE; 6126f24a144SLois Curfman McInnes snes->vwork = 0; 6136f24a144SLois Curfman McInnes snes->nwork = 0; 614758f92a0SBarry Smith snes->conv_hist_len = 0; 615758f92a0SBarry Smith snes->conv_hist_max = 0; 616758f92a0SBarry Smith snes->conv_hist = PETSC_NULL; 617758f92a0SBarry Smith snes->conv_hist_its = PETSC_NULL; 618758f92a0SBarry Smith snes->conv_hist_reset = PETSC_TRUE; 619184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 6209b94acceSBarry Smith 6219b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 622b0a32e0cSBarry Smith ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 62352e6d16bSBarry Smith ierr = PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx));CHKERRQ(ierr); 6249b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6259b94acceSBarry Smith kctx->version = 2; 6269b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6279b94acceSBarry Smith this was too large for some test cases */ 6289b94acceSBarry Smith kctx->rtol_last = 0; 6299b94acceSBarry Smith kctx->rtol_max = .9; 6309b94acceSBarry Smith kctx->gamma = 1.0; 6319b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6329b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6339b94acceSBarry Smith kctx->threshold = .1; 6349b94acceSBarry Smith kctx->lresid_last = 0; 6359b94acceSBarry Smith kctx->norm_last = 0; 6369b94acceSBarry Smith 63794b7f48cSBarry Smith ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr); 63852e6d16bSBarry Smith ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr); 6395334005bSBarry Smith 6409b94acceSBarry Smith *outsnes = snes; 64100036973SBarry Smith ierr = PetscPublishAll(snes);CHKERRQ(ierr); 6423a40ed3dSBarry Smith PetscFunctionReturn(0); 6439b94acceSBarry Smith } 6449b94acceSBarry Smith 6454a2ae208SSatish Balay #undef __FUNCT__ 6464a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction" 6479b94acceSBarry Smith /*@C 6489b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6499b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6509b94acceSBarry Smith equations. 6519b94acceSBarry Smith 652fee21e36SBarry Smith Collective on SNES 653fee21e36SBarry Smith 654c7afd0dbSLois Curfman McInnes Input Parameters: 655c7afd0dbSLois Curfman McInnes + snes - the SNES context 656c7afd0dbSLois Curfman McInnes . func - function evaluation routine 657c7afd0dbSLois Curfman McInnes . r - vector to store function value 658c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 659c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6609b94acceSBarry Smith 661c7afd0dbSLois Curfman McInnes Calling sequence of func: 6628d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 663c7afd0dbSLois Curfman McInnes 664313e4042SLois Curfman McInnes . f - function vector 665c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 6669b94acceSBarry Smith 6679b94acceSBarry Smith Notes: 6689b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6699b94acceSBarry Smith $ f'(x) x = -f(x), 670c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 6719b94acceSBarry Smith 67236851e7fSLois Curfman McInnes Level: beginner 67336851e7fSLois Curfman McInnes 6749b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6759b94acceSBarry Smith 676a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6779b94acceSBarry Smith @*/ 6786849ba73SBarry Smith PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 6799b94acceSBarry Smith { 6803a40ed3dSBarry Smith PetscFunctionBegin; 6814482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 6824482741eSBarry Smith PetscValidHeaderSpecific(r,VEC_COOKIE,2); 683c9780b6fSBarry Smith PetscCheckSameComm(snes,1,r,2); 684184914b5SBarry Smith 6859b94acceSBarry Smith snes->computefunction = func; 6869b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6879b94acceSBarry Smith snes->funP = ctx; 6883a40ed3dSBarry Smith PetscFunctionReturn(0); 6899b94acceSBarry Smith } 6909b94acceSBarry Smith 6913ab0aad5SBarry Smith /* --------------------------------------------------------------- */ 6923ab0aad5SBarry Smith #undef __FUNCT__ 6933ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs" 6943ab0aad5SBarry Smith /*@C 6953ab0aad5SBarry Smith SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set 6963ab0aad5SBarry Smith it assumes a zero right hand side. 6973ab0aad5SBarry Smith 6983ab0aad5SBarry Smith Collective on SNES 6993ab0aad5SBarry Smith 7003ab0aad5SBarry Smith Input Parameters: 7013ab0aad5SBarry Smith + snes - the SNES context 7023ab0aad5SBarry Smith - rhs - the right hand side vector or PETSC_NULL for a zero right hand side 7033ab0aad5SBarry Smith 7043ab0aad5SBarry Smith Level: intermediate 7053ab0aad5SBarry Smith 7063ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side 7073ab0aad5SBarry Smith 7083ab0aad5SBarry Smith .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 7093ab0aad5SBarry Smith @*/ 710dfbe8321SBarry Smith PetscErrorCode SNESSetRhs(SNES snes,Vec rhs) 7113ab0aad5SBarry Smith { 712dfbe8321SBarry Smith PetscErrorCode ierr; 7133ab0aad5SBarry Smith 7143ab0aad5SBarry Smith PetscFunctionBegin; 7153ab0aad5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7163ab0aad5SBarry Smith if (rhs) { 7173ab0aad5SBarry Smith PetscValidHeaderSpecific(rhs,VEC_COOKIE,2); 7183ab0aad5SBarry Smith PetscCheckSameComm(snes,1,rhs,2); 7193ab0aad5SBarry Smith ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr); 7203ab0aad5SBarry Smith } 7213ab0aad5SBarry Smith if (snes->afine) { 7223ab0aad5SBarry Smith ierr = VecDestroy(snes->afine);CHKERRQ(ierr); 7233ab0aad5SBarry Smith } 7243ab0aad5SBarry Smith snes->afine = rhs; 7253ab0aad5SBarry Smith PetscFunctionReturn(0); 7263ab0aad5SBarry Smith } 7273ab0aad5SBarry Smith 7284a2ae208SSatish Balay #undef __FUNCT__ 7294a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction" 7309b94acceSBarry Smith /*@ 73136851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7329b94acceSBarry Smith SNESSetFunction(). 7339b94acceSBarry Smith 734c7afd0dbSLois Curfman McInnes Collective on SNES 735c7afd0dbSLois Curfman McInnes 7369b94acceSBarry Smith Input Parameters: 737c7afd0dbSLois Curfman McInnes + snes - the SNES context 738c7afd0dbSLois Curfman McInnes - x - input vector 7399b94acceSBarry Smith 7409b94acceSBarry Smith Output Parameter: 7413638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7429b94acceSBarry Smith 7431bffabb2SLois Curfman McInnes Notes: 74436851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 74536851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 74636851e7fSLois Curfman McInnes themselves. 74736851e7fSLois Curfman McInnes 74836851e7fSLois Curfman McInnes Level: developer 74936851e7fSLois Curfman McInnes 7509b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 7519b94acceSBarry Smith 752a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 7539b94acceSBarry Smith @*/ 754dfbe8321SBarry Smith PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y) 7559b94acceSBarry Smith { 756dfbe8321SBarry Smith PetscErrorCode ierr; 7579b94acceSBarry Smith 7583a40ed3dSBarry Smith PetscFunctionBegin; 7594482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7604482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 7614482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,3); 762c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 763c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,3); 764184914b5SBarry Smith 765d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 766d64ed03dSBarry Smith PetscStackPush("SNES user function"); 7679b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 768d64ed03dSBarry Smith PetscStackPop; 7693ab0aad5SBarry Smith if (snes->afine) { 7703ab0aad5SBarry Smith PetscScalar mone = -1.0; 7713ab0aad5SBarry Smith ierr = VecAXPY(&mone,snes->afine,y);CHKERRQ(ierr); 7723ab0aad5SBarry Smith } 773ae3c334cSLois Curfman McInnes snes->nfuncs++; 774d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 7753a40ed3dSBarry Smith PetscFunctionReturn(0); 7769b94acceSBarry Smith } 7779b94acceSBarry Smith 7784a2ae208SSatish Balay #undef __FUNCT__ 7794a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian" 78062fef451SLois Curfman McInnes /*@ 78162fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 78262fef451SLois Curfman McInnes set with SNESSetJacobian(). 78362fef451SLois Curfman McInnes 784c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 785c7afd0dbSLois Curfman McInnes 78662fef451SLois Curfman McInnes Input Parameters: 787c7afd0dbSLois Curfman McInnes + snes - the SNES context 788c7afd0dbSLois Curfman McInnes - x - input vector 78962fef451SLois Curfman McInnes 79062fef451SLois Curfman McInnes Output Parameters: 791c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 79262fef451SLois Curfman McInnes . B - optional preconditioning matrix 793c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 794fee21e36SBarry Smith 79562fef451SLois Curfman McInnes Notes: 79662fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 79762fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 79862fef451SLois Curfman McInnes 79994b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 800dc5a77f8SLois Curfman McInnes flag parameter. 80162fef451SLois Curfman McInnes 80236851e7fSLois Curfman McInnes Level: developer 80336851e7fSLois Curfman McInnes 80462fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 80562fef451SLois Curfman McInnes 80694b7f48cSBarry Smith .seealso: SNESSetJacobian(), KSPSetOperators() 80762fef451SLois Curfman McInnes @*/ 808dfbe8321SBarry Smith PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8099b94acceSBarry Smith { 810dfbe8321SBarry Smith PetscErrorCode ierr; 8113a40ed3dSBarry Smith 8123a40ed3dSBarry Smith PetscFunctionBegin; 8134482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8144482741eSBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE,2); 8154482741eSBarry Smith PetscValidPointer(flg,5); 816c9780b6fSBarry Smith PetscCheckSameComm(snes,1,X,2); 8173a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 818d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 819c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 820d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 8219b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 822d64ed03dSBarry Smith PetscStackPop; 823d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 8246d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 8254482741eSBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 8264482741eSBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 8273a40ed3dSBarry Smith PetscFunctionReturn(0); 8289b94acceSBarry Smith } 8299b94acceSBarry Smith 8304a2ae208SSatish Balay #undef __FUNCT__ 8314a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian" 8329b94acceSBarry Smith /*@C 8339b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 834044dda88SLois Curfman McInnes location to store the matrix. 8359b94acceSBarry Smith 836c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 837c7afd0dbSLois Curfman McInnes 8389b94acceSBarry Smith Input Parameters: 839c7afd0dbSLois Curfman McInnes + snes - the SNES context 8409b94acceSBarry Smith . A - Jacobian matrix 8419b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 8429b94acceSBarry Smith . func - Jacobian evaluation routine 843c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 8442cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 8459b94acceSBarry Smith 8469b94acceSBarry Smith Calling sequence of func: 8478d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 8489b94acceSBarry Smith 849c7afd0dbSLois Curfman McInnes + x - input vector 8509b94acceSBarry Smith . A - Jacobian matrix 8519b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 852ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 85394b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 854c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 8559b94acceSBarry Smith 8569b94acceSBarry Smith Notes: 85794b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 8582cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 859ac21db08SLois Curfman McInnes 860ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 8619b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 8629b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 8639b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 8649b94acceSBarry Smith throughout the global iterations. 8659b94acceSBarry Smith 86636851e7fSLois Curfman McInnes Level: beginner 86736851e7fSLois Curfman McInnes 8689b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 8699b94acceSBarry Smith 87013f74950SBarry Smith .seealso: KSPSetOperators(), SNESSetFunction(), , MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor() 8719b94acceSBarry Smith @*/ 8726849ba73SBarry Smith PetscErrorCode SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 8739b94acceSBarry Smith { 874dfbe8321SBarry Smith PetscErrorCode ierr; 8753a7fca6bSBarry Smith 8763a40ed3dSBarry Smith PetscFunctionBegin; 8774482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8784482741eSBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 8794482741eSBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 880c9780b6fSBarry Smith if (A) PetscCheckSameComm(snes,1,A,2); 881c9780b6fSBarry Smith if (B) PetscCheckSameComm(snes,1,B,2); 8823a7fca6bSBarry Smith if (func) snes->computejacobian = func; 8833a7fca6bSBarry Smith if (ctx) snes->jacP = ctx; 8843a7fca6bSBarry Smith if (A) { 8853a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 8869b94acceSBarry Smith snes->jacobian = A; 8873a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 8883a7fca6bSBarry Smith } 8893a7fca6bSBarry Smith if (B) { 8903a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 8919b94acceSBarry Smith snes->jacobian_pre = B; 8923a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 8933a7fca6bSBarry Smith } 89469a612a9SBarry Smith ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 8953a40ed3dSBarry Smith PetscFunctionReturn(0); 8969b94acceSBarry Smith } 89762fef451SLois Curfman McInnes 8984a2ae208SSatish Balay #undef __FUNCT__ 8994a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian" 900c2aafc4cSSatish Balay /*@C 901b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 902b4fd4287SBarry Smith provided context for evaluating the Jacobian. 903b4fd4287SBarry Smith 904c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 905c7afd0dbSLois Curfman McInnes 906b4fd4287SBarry Smith Input Parameter: 907b4fd4287SBarry Smith . snes - the nonlinear solver context 908b4fd4287SBarry Smith 909b4fd4287SBarry Smith Output Parameters: 910c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 911b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 91200036973SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 91300036973SBarry Smith - func - location to put Jacobian function (or PETSC_NULL) 914fee21e36SBarry Smith 91536851e7fSLois Curfman McInnes Level: advanced 91636851e7fSLois Curfman McInnes 917b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 918b4fd4287SBarry Smith @*/ 9196849ba73SBarry Smith PetscErrorCode SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 920b4fd4287SBarry Smith { 9213a40ed3dSBarry Smith PetscFunctionBegin; 9224482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 923b4fd4287SBarry Smith if (A) *A = snes->jacobian; 924b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 925b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 92600036973SBarry Smith if (func) *func = snes->computejacobian; 9273a40ed3dSBarry Smith PetscFunctionReturn(0); 928b4fd4287SBarry Smith } 929b4fd4287SBarry Smith 9309b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 931dfbe8321SBarry Smith EXTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 9329b94acceSBarry Smith 9334a2ae208SSatish Balay #undef __FUNCT__ 9344a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp" 9359b94acceSBarry Smith /*@ 9369b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 937272ac6f2SLois Curfman McInnes of a nonlinear solver. 9389b94acceSBarry Smith 939fee21e36SBarry Smith Collective on SNES 940fee21e36SBarry Smith 941c7afd0dbSLois Curfman McInnes Input Parameters: 942c7afd0dbSLois Curfman McInnes + snes - the SNES context 943c7afd0dbSLois Curfman McInnes - x - the solution vector 944c7afd0dbSLois Curfman McInnes 945272ac6f2SLois Curfman McInnes Notes: 946272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 947272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 948272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 949272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 950272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 951272ac6f2SLois Curfman McInnes 95236851e7fSLois Curfman McInnes Level: advanced 95336851e7fSLois Curfman McInnes 9549b94acceSBarry Smith .keywords: SNES, nonlinear, setup 9559b94acceSBarry Smith 9569b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 9579b94acceSBarry Smith @*/ 958dfbe8321SBarry Smith PetscErrorCode SNESSetUp(SNES snes,Vec x) 9599b94acceSBarry Smith { 960dfbe8321SBarry Smith PetscErrorCode ierr; 9614b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 9623a40ed3dSBarry Smith 9633a40ed3dSBarry Smith PetscFunctionBegin; 9644482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 9654482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 966c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 9678ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 9689b94acceSBarry Smith 969b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 970c1f60f51SBarry Smith /* 971c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 972dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 973c1f60f51SBarry Smith */ 974c1f60f51SBarry Smith if (flg) { 975c1f60f51SBarry Smith Mat J; 9765a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 9775a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 9783a7fca6bSBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n"); 9793a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 9803a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 981c1f60f51SBarry Smith } 98245fc7adcSBarry Smith 983b4014bb2SMatthew Knepley #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) 98445fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 98545fc7adcSBarry Smith if (flg) { 98645fc7adcSBarry Smith Mat J; 98745fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 98845fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 98945fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 99045fc7adcSBarry Smith } 99132a4b47aSMatthew Knepley #endif 99245fc7adcSBarry Smith 993b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 994c1f60f51SBarry Smith /* 995dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 996c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 997c1f60f51SBarry Smith */ 99831615d04SLois Curfman McInnes if (flg) { 999272ac6f2SLois Curfman McInnes Mat J; 1000b5d62d44SBarry Smith KSP ksp; 100194b7f48cSBarry Smith PC pc; 1002f3ef73ceSBarry Smith 10035a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10043a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 100593b98531SBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"); 10063a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 10073a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 10083a7fca6bSBarry Smith 1009f3ef73ceSBarry Smith /* force no preconditioner */ 101094b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1011b5d62d44SBarry Smith ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1012a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1013a9815358SBarry Smith if (!flg) { 1014f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1015272ac6f2SLois Curfman McInnes } 1016a9815358SBarry Smith } 1017f3ef73ceSBarry Smith 101829bbc08cSBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 101929bbc08cSBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 102029bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1021a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 102229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1023a8c6a408SBarry Smith } 1024a703fe33SLois Curfman McInnes 1025a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 10264b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 10276831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 102894b7f48cSBarry Smith KSP ksp; 102994b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1030186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1031a703fe33SLois Curfman McInnes } 10324b27c08aSLois Curfman McInnes 1033a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 103482bf6240SBarry Smith snes->setupcalled = 1; 10353a40ed3dSBarry Smith PetscFunctionReturn(0); 10369b94acceSBarry Smith } 10379b94acceSBarry Smith 10384a2ae208SSatish Balay #undef __FUNCT__ 10394a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 10409b94acceSBarry Smith /*@C 10419b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 10429b94acceSBarry Smith with SNESCreate(). 10439b94acceSBarry Smith 1044c7afd0dbSLois Curfman McInnes Collective on SNES 1045c7afd0dbSLois Curfman McInnes 10469b94acceSBarry Smith Input Parameter: 10479b94acceSBarry Smith . snes - the SNES context 10489b94acceSBarry Smith 104936851e7fSLois Curfman McInnes Level: beginner 105036851e7fSLois Curfman McInnes 10519b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 10529b94acceSBarry Smith 105363a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 10549b94acceSBarry Smith @*/ 1055dfbe8321SBarry Smith PetscErrorCode SNESDestroy(SNES snes) 10569b94acceSBarry Smith { 105777431f27SBarry Smith PetscInt i; 10586849ba73SBarry Smith PetscErrorCode ierr; 10593a40ed3dSBarry Smith 10603a40ed3dSBarry Smith PetscFunctionBegin; 10614482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 10623a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1063d4bb536fSBarry Smith 1064be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 10650f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1066be0abb6dSBarry Smith 1067e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1068606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 10693a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 10703a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 10713ab0aad5SBarry Smith if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 107294b7f48cSBarry Smith ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1073522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1074b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1075b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1076b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1077b8a78c4aSBarry Smith } 1078b8a78c4aSBarry Smith } 1079*a79aaaedSSatish Balay ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 10803a40ed3dSBarry Smith PetscFunctionReturn(0); 10819b94acceSBarry Smith } 10829b94acceSBarry Smith 10839b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 10849b94acceSBarry Smith 10854a2ae208SSatish Balay #undef __FUNCT__ 10864a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 10879b94acceSBarry Smith /*@ 1088d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 10899b94acceSBarry Smith 1090c7afd0dbSLois Curfman McInnes Collective on SNES 1091c7afd0dbSLois Curfman McInnes 10929b94acceSBarry Smith Input Parameters: 1093c7afd0dbSLois Curfman McInnes + snes - the SNES context 109470441072SBarry Smith . abstol - absolute convergence tolerance 109533174efeSLois Curfman McInnes . rtol - relative convergence tolerance 109633174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 109733174efeSLois Curfman McInnes of the change in the solution between steps 109833174efeSLois Curfman McInnes . maxit - maximum number of iterations 1099c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1100fee21e36SBarry Smith 110133174efeSLois Curfman McInnes Options Database Keys: 110270441072SBarry Smith + -snes_atol <abstol> - Sets abstol 1103c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1104c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1105c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1106c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 11079b94acceSBarry Smith 1108d7a720efSLois Curfman McInnes Notes: 11099b94acceSBarry Smith The default maximum number of iterations is 50. 11109b94acceSBarry Smith The default maximum number of function evaluations is 1000. 11119b94acceSBarry Smith 111236851e7fSLois Curfman McInnes Level: intermediate 111336851e7fSLois Curfman McInnes 111433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 11159b94acceSBarry Smith 11162492ecdbSBarry Smith .seealso: SNESSetTrustRegionTolerance() 11179b94acceSBarry Smith @*/ 111870441072SBarry Smith PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 11199b94acceSBarry Smith { 11203a40ed3dSBarry Smith PetscFunctionBegin; 11214482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 112270441072SBarry Smith if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1123d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1124d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1125d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1126d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11273a40ed3dSBarry Smith PetscFunctionReturn(0); 11289b94acceSBarry Smith } 11299b94acceSBarry Smith 11304a2ae208SSatish Balay #undef __FUNCT__ 11314a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11329b94acceSBarry Smith /*@ 113333174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 113433174efeSLois Curfman McInnes 1135c7afd0dbSLois Curfman McInnes Not Collective 1136c7afd0dbSLois Curfman McInnes 113733174efeSLois Curfman McInnes Input Parameters: 1138c7afd0dbSLois Curfman McInnes + snes - the SNES context 113970441072SBarry Smith . abstol - absolute convergence tolerance 114033174efeSLois Curfman McInnes . rtol - relative convergence tolerance 114133174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 114233174efeSLois Curfman McInnes of the change in the solution between steps 114333174efeSLois Curfman McInnes . maxit - maximum number of iterations 1144c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1145fee21e36SBarry Smith 114633174efeSLois Curfman McInnes Notes: 114733174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 114833174efeSLois Curfman McInnes 114936851e7fSLois Curfman McInnes Level: intermediate 115036851e7fSLois Curfman McInnes 115133174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 115233174efeSLois Curfman McInnes 115333174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 115433174efeSLois Curfman McInnes @*/ 115570441072SBarry Smith PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 115633174efeSLois Curfman McInnes { 11573a40ed3dSBarry Smith PetscFunctionBegin; 11584482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 115970441072SBarry Smith if (abstol) *abstol = snes->abstol; 116033174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 116133174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 116233174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 116333174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 11643a40ed3dSBarry Smith PetscFunctionReturn(0); 116533174efeSLois Curfman McInnes } 116633174efeSLois Curfman McInnes 11674a2ae208SSatish Balay #undef __FUNCT__ 11684a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 116933174efeSLois Curfman McInnes /*@ 11709b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 11719b94acceSBarry Smith 1172fee21e36SBarry Smith Collective on SNES 1173fee21e36SBarry Smith 1174c7afd0dbSLois Curfman McInnes Input Parameters: 1175c7afd0dbSLois Curfman McInnes + snes - the SNES context 1176c7afd0dbSLois Curfman McInnes - tol - tolerance 1177c7afd0dbSLois Curfman McInnes 11789b94acceSBarry Smith Options Database Key: 1179c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 11809b94acceSBarry Smith 118136851e7fSLois Curfman McInnes Level: intermediate 118236851e7fSLois Curfman McInnes 11839b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 11849b94acceSBarry Smith 11852492ecdbSBarry Smith .seealso: SNESSetTolerances() 11869b94acceSBarry Smith @*/ 1187dfbe8321SBarry Smith PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 11889b94acceSBarry Smith { 11893a40ed3dSBarry Smith PetscFunctionBegin; 11904482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 11919b94acceSBarry Smith snes->deltatol = tol; 11923a40ed3dSBarry Smith PetscFunctionReturn(0); 11939b94acceSBarry Smith } 11949b94acceSBarry Smith 1195df9fa365SBarry Smith /* 1196df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1197df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1198df9fa365SBarry Smith macros instead of functions 1199df9fa365SBarry Smith */ 12004a2ae208SSatish Balay #undef __FUNCT__ 12014a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 1202a7cc72afSBarry Smith PetscErrorCode SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1203ce1608b8SBarry Smith { 1204dfbe8321SBarry Smith PetscErrorCode ierr; 1205ce1608b8SBarry Smith 1206ce1608b8SBarry Smith PetscFunctionBegin; 12074482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1208ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1209ce1608b8SBarry Smith PetscFunctionReturn(0); 1210ce1608b8SBarry Smith } 1211ce1608b8SBarry Smith 12124a2ae208SSatish Balay #undef __FUNCT__ 12134a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 1214dfbe8321SBarry Smith PetscErrorCode SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1215df9fa365SBarry Smith { 1216dfbe8321SBarry Smith PetscErrorCode ierr; 1217df9fa365SBarry Smith 1218df9fa365SBarry Smith PetscFunctionBegin; 1219df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1220df9fa365SBarry Smith PetscFunctionReturn(0); 1221df9fa365SBarry Smith } 1222df9fa365SBarry Smith 12234a2ae208SSatish Balay #undef __FUNCT__ 12244a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 1225dfbe8321SBarry Smith PetscErrorCode SNESLGMonitorDestroy(PetscDrawLG draw) 1226df9fa365SBarry Smith { 1227dfbe8321SBarry Smith PetscErrorCode ierr; 1228df9fa365SBarry Smith 1229df9fa365SBarry Smith PetscFunctionBegin; 1230df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1231df9fa365SBarry Smith PetscFunctionReturn(0); 1232df9fa365SBarry Smith } 1233df9fa365SBarry Smith 12349b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12359b94acceSBarry Smith 12364a2ae208SSatish Balay #undef __FUNCT__ 12374a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12389b94acceSBarry Smith /*@C 12395cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12409b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12419b94acceSBarry Smith progress. 12429b94acceSBarry Smith 1243fee21e36SBarry Smith Collective on SNES 1244fee21e36SBarry Smith 1245c7afd0dbSLois Curfman McInnes Input Parameters: 1246c7afd0dbSLois Curfman McInnes + snes - the SNES context 1247c7afd0dbSLois Curfman McInnes . func - monitoring routine 1248b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1249b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desitre) 1250b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1251b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 12529b94acceSBarry Smith 1253c7afd0dbSLois Curfman McInnes Calling sequence of func: 1254a7cc72afSBarry Smith $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1255c7afd0dbSLois Curfman McInnes 1256c7afd0dbSLois Curfman McInnes + snes - the SNES context 1257c7afd0dbSLois Curfman McInnes . its - iteration number 1258c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 125940a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 12609b94acceSBarry Smith 12619665c990SLois Curfman McInnes Options Database Keys: 1262c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1263c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1264c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1265c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1266c7afd0dbSLois Curfman McInnes been hardwired into a code by 1267c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1268c7afd0dbSLois Curfman McInnes does not cancel those set via 1269c7afd0dbSLois Curfman McInnes the options database. 12709665c990SLois Curfman McInnes 1271639f9d9dSBarry Smith Notes: 12726bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 12736bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 12746bc08f3fSLois Curfman McInnes order in which they were set. 1275639f9d9dSBarry Smith 127636851e7fSLois Curfman McInnes Level: intermediate 127736851e7fSLois Curfman McInnes 12789b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 12799b94acceSBarry Smith 12805cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 12819b94acceSBarry Smith @*/ 128277431f27SBarry Smith PetscErrorCode SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 12839b94acceSBarry Smith { 12843a40ed3dSBarry Smith PetscFunctionBegin; 12854482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1286639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 128729bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1288639f9d9dSBarry Smith } 1289639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1290b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1291639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 12923a40ed3dSBarry Smith PetscFunctionReturn(0); 12939b94acceSBarry Smith } 12949b94acceSBarry Smith 12954a2ae208SSatish Balay #undef __FUNCT__ 12964a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 12975cd90555SBarry Smith /*@C 12985cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 12995cd90555SBarry Smith 1300c7afd0dbSLois Curfman McInnes Collective on SNES 1301c7afd0dbSLois Curfman McInnes 13025cd90555SBarry Smith Input Parameters: 13035cd90555SBarry Smith . snes - the SNES context 13045cd90555SBarry Smith 13051a480d89SAdministrator Options Database Key: 1306c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1307c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1308c7afd0dbSLois Curfman McInnes set via the options database 13095cd90555SBarry Smith 13105cd90555SBarry Smith Notes: 13115cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 13125cd90555SBarry Smith 131336851e7fSLois Curfman McInnes Level: intermediate 131436851e7fSLois Curfman McInnes 13155cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 13165cd90555SBarry Smith 13175cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 13185cd90555SBarry Smith @*/ 1319dfbe8321SBarry Smith PetscErrorCode SNESClearMonitor(SNES snes) 13205cd90555SBarry Smith { 13215cd90555SBarry Smith PetscFunctionBegin; 13224482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13235cd90555SBarry Smith snes->numbermonitors = 0; 13245cd90555SBarry Smith PetscFunctionReturn(0); 13255cd90555SBarry Smith } 13265cd90555SBarry Smith 13274a2ae208SSatish Balay #undef __FUNCT__ 13284a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13299b94acceSBarry Smith /*@C 13309b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13319b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13329b94acceSBarry Smith 1333fee21e36SBarry Smith Collective on SNES 1334fee21e36SBarry Smith 1335c7afd0dbSLois Curfman McInnes Input Parameters: 1336c7afd0dbSLois Curfman McInnes + snes - the SNES context 1337c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1338c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1339c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 13409b94acceSBarry Smith 1341c7afd0dbSLois Curfman McInnes Calling sequence of func: 13426849ba73SBarry Smith $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1343c7afd0dbSLois Curfman McInnes 1344c7afd0dbSLois Curfman McInnes + snes - the SNES context 1345c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1346184914b5SBarry Smith . reason - reason for convergence/divergence 1347c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 13484b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 13494b27c08aSLois Curfman McInnes - f - 2-norm of function 13509b94acceSBarry Smith 135136851e7fSLois Curfman McInnes Level: advanced 135236851e7fSLois Curfman McInnes 13539b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13549b94acceSBarry Smith 13554b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 13569b94acceSBarry Smith @*/ 13576849ba73SBarry Smith PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 13589b94acceSBarry Smith { 13593a40ed3dSBarry Smith PetscFunctionBegin; 13604482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13619b94acceSBarry Smith (snes)->converged = func; 13629b94acceSBarry Smith (snes)->cnvP = cctx; 13633a40ed3dSBarry Smith PetscFunctionReturn(0); 13649b94acceSBarry Smith } 13659b94acceSBarry Smith 13664a2ae208SSatish Balay #undef __FUNCT__ 13674a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 1368184914b5SBarry Smith /*@C 1369184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1370184914b5SBarry Smith 1371184914b5SBarry Smith Not Collective 1372184914b5SBarry Smith 1373184914b5SBarry Smith Input Parameter: 1374184914b5SBarry Smith . snes - the SNES context 1375184914b5SBarry Smith 1376184914b5SBarry Smith Output Parameter: 1377e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1378184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1379184914b5SBarry Smith 1380184914b5SBarry Smith Level: intermediate 1381184914b5SBarry Smith 1382184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1383184914b5SBarry Smith 1384184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1385184914b5SBarry Smith 13864b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1387184914b5SBarry Smith @*/ 1388dfbe8321SBarry Smith PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1389184914b5SBarry Smith { 1390184914b5SBarry Smith PetscFunctionBegin; 13914482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13924482741eSBarry Smith PetscValidPointer(reason,2); 1393184914b5SBarry Smith *reason = snes->reason; 1394184914b5SBarry Smith PetscFunctionReturn(0); 1395184914b5SBarry Smith } 1396184914b5SBarry Smith 13974a2ae208SSatish Balay #undef __FUNCT__ 13984a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1399c9005455SLois Curfman McInnes /*@ 1400c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1401c9005455SLois Curfman McInnes 1402fee21e36SBarry Smith Collective on SNES 1403fee21e36SBarry Smith 1404c7afd0dbSLois Curfman McInnes Input Parameters: 1405c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1406c7afd0dbSLois Curfman McInnes . a - array to hold history 1407cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1408758f92a0SBarry Smith . na - size of a and its 140964731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1410758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1411c7afd0dbSLois Curfman McInnes 1412c9005455SLois Curfman McInnes Notes: 14134b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1414c9005455SLois Curfman McInnes at each step. 1415c9005455SLois Curfman McInnes 1416c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1417c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1418c9005455SLois Curfman McInnes during the section of code that is being timed. 1419c9005455SLois Curfman McInnes 142036851e7fSLois Curfman McInnes Level: intermediate 142136851e7fSLois Curfman McInnes 1422c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1423758f92a0SBarry Smith 142408405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1425758f92a0SBarry Smith 1426c9005455SLois Curfman McInnes @*/ 142777431f27SBarry Smith PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1428c9005455SLois Curfman McInnes { 14293a40ed3dSBarry Smith PetscFunctionBegin; 14304482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14314482741eSBarry Smith if (na) PetscValidScalarPointer(a,2); 1432c9005455SLois Curfman McInnes snes->conv_hist = a; 1433758f92a0SBarry Smith snes->conv_hist_its = its; 1434758f92a0SBarry Smith snes->conv_hist_max = na; 1435758f92a0SBarry Smith snes->conv_hist_reset = reset; 1436758f92a0SBarry Smith PetscFunctionReturn(0); 1437758f92a0SBarry Smith } 1438758f92a0SBarry Smith 14394a2ae208SSatish Balay #undef __FUNCT__ 14404a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 14410c4c9dddSBarry Smith /*@C 1442758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1443758f92a0SBarry Smith 1444758f92a0SBarry Smith Collective on SNES 1445758f92a0SBarry Smith 1446758f92a0SBarry Smith Input Parameter: 1447758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1448758f92a0SBarry Smith 1449758f92a0SBarry Smith Output Parameters: 1450758f92a0SBarry Smith . a - array to hold history 1451758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1452758f92a0SBarry Smith negative if not converged) for each solve. 1453758f92a0SBarry Smith - na - size of a and its 1454758f92a0SBarry Smith 1455758f92a0SBarry Smith Notes: 1456758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1457758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1458758f92a0SBarry Smith 1459758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1460758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1461758f92a0SBarry Smith during the section of code that is being timed. 1462758f92a0SBarry Smith 1463758f92a0SBarry Smith Level: intermediate 1464758f92a0SBarry Smith 1465758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1466758f92a0SBarry Smith 1467758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1468758f92a0SBarry Smith 1469758f92a0SBarry Smith @*/ 147077431f27SBarry Smith PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1471758f92a0SBarry Smith { 1472758f92a0SBarry Smith PetscFunctionBegin; 14734482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1474758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1475758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1476758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 14773a40ed3dSBarry Smith PetscFunctionReturn(0); 1478c9005455SLois Curfman McInnes } 1479c9005455SLois Curfman McInnes 1480e74ef692SMatthew Knepley #undef __FUNCT__ 1481e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC" 1482ac226902SBarry Smith /*@C 148376b2cf59SMatthew Knepley SNESSetRhsBC - Sets the function which applies boundary conditions 148476b2cf59SMatthew Knepley to the Rhs of each system. 148576b2cf59SMatthew Knepley 148676b2cf59SMatthew Knepley Collective on SNES 148776b2cf59SMatthew Knepley 148876b2cf59SMatthew Knepley Input Parameters: 148976b2cf59SMatthew Knepley . snes - The nonlinear solver context 149076b2cf59SMatthew Knepley . func - The function 149176b2cf59SMatthew Knepley 149276b2cf59SMatthew Knepley Calling sequence of func: 149376b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx); 149476b2cf59SMatthew Knepley 149576b2cf59SMatthew Knepley . rhs - The current rhs vector 149676b2cf59SMatthew Knepley . ctx - The user-context 149776b2cf59SMatthew Knepley 149876b2cf59SMatthew Knepley Level: intermediate 149976b2cf59SMatthew Knepley 150076b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 150176b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate() 150276b2cf59SMatthew Knepley @*/ 15036849ba73SBarry Smith PetscErrorCode SNESSetRhsBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 150476b2cf59SMatthew Knepley { 150576b2cf59SMatthew Knepley PetscFunctionBegin; 15064482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 150776b2cf59SMatthew Knepley snes->applyrhsbc = func; 150876b2cf59SMatthew Knepley PetscFunctionReturn(0); 150976b2cf59SMatthew Knepley } 151076b2cf59SMatthew Knepley 1511e74ef692SMatthew Knepley #undef __FUNCT__ 1512e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC" 151376b2cf59SMatthew Knepley /*@ 151476b2cf59SMatthew Knepley SNESDefaultRhsBC - The default boundary condition function which does nothing. 151576b2cf59SMatthew Knepley 151676b2cf59SMatthew Knepley Not collective 151776b2cf59SMatthew Knepley 151876b2cf59SMatthew Knepley Input Parameters: 151976b2cf59SMatthew Knepley . snes - The nonlinear solver context 152076b2cf59SMatthew Knepley . rhs - The Rhs 152176b2cf59SMatthew Knepley . ctx - The user-context 152276b2cf59SMatthew Knepley 152376b2cf59SMatthew Knepley Level: beginner 152476b2cf59SMatthew Knepley 152576b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 152676b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate() 152776b2cf59SMatthew Knepley @*/ 1528dfbe8321SBarry Smith PetscErrorCode SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx) 152976b2cf59SMatthew Knepley { 153076b2cf59SMatthew Knepley PetscFunctionBegin; 153176b2cf59SMatthew Knepley PetscFunctionReturn(0); 153276b2cf59SMatthew Knepley } 153376b2cf59SMatthew Knepley 1534e74ef692SMatthew Knepley #undef __FUNCT__ 1535e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC" 1536ac226902SBarry Smith /*@C 153776b2cf59SMatthew Knepley SNESSetSolutionBC - Sets the function which applies boundary conditions 153876b2cf59SMatthew Knepley to the solution of each system. 153976b2cf59SMatthew Knepley 154076b2cf59SMatthew Knepley Collective on SNES 154176b2cf59SMatthew Knepley 154276b2cf59SMatthew Knepley Input Parameters: 154376b2cf59SMatthew Knepley . snes - The nonlinear solver context 154476b2cf59SMatthew Knepley . func - The function 154576b2cf59SMatthew Knepley 154676b2cf59SMatthew Knepley Calling sequence of func: 15479d06fb6bSMatthew Knepley . func (SNES snes, Vec rsol, void *ctx); 154876b2cf59SMatthew Knepley 154976b2cf59SMatthew Knepley . sol - The current solution vector 155076b2cf59SMatthew Knepley . ctx - The user-context 155176b2cf59SMatthew Knepley 155276b2cf59SMatthew Knepley Level: intermediate 155376b2cf59SMatthew Knepley 155476b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 155576b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate() 155676b2cf59SMatthew Knepley @*/ 15576849ba73SBarry Smith PetscErrorCode SNESSetSolutionBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 155876b2cf59SMatthew Knepley { 155976b2cf59SMatthew Knepley PetscFunctionBegin; 15604482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 156176b2cf59SMatthew Knepley snes->applysolbc = func; 156276b2cf59SMatthew Knepley PetscFunctionReturn(0); 156376b2cf59SMatthew Knepley } 156476b2cf59SMatthew Knepley 1565e74ef692SMatthew Knepley #undef __FUNCT__ 1566e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC" 156776b2cf59SMatthew Knepley /*@ 156876b2cf59SMatthew Knepley SNESDefaultSolutionBC - The default boundary condition function which does nothing. 156976b2cf59SMatthew Knepley 157076b2cf59SMatthew Knepley Not collective 157176b2cf59SMatthew Knepley 157276b2cf59SMatthew Knepley Input Parameters: 157376b2cf59SMatthew Knepley . snes - The nonlinear solver context 157476b2cf59SMatthew Knepley . sol - The solution 157576b2cf59SMatthew Knepley . ctx - The user-context 157676b2cf59SMatthew Knepley 157776b2cf59SMatthew Knepley Level: beginner 157876b2cf59SMatthew Knepley 157976b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 158076b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate() 158176b2cf59SMatthew Knepley @*/ 1582dfbe8321SBarry Smith PetscErrorCode SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx) 158376b2cf59SMatthew Knepley { 158476b2cf59SMatthew Knepley PetscFunctionBegin; 158576b2cf59SMatthew Knepley PetscFunctionReturn(0); 158676b2cf59SMatthew Knepley } 158776b2cf59SMatthew Knepley 1588e74ef692SMatthew Knepley #undef __FUNCT__ 1589e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 1590ac226902SBarry Smith /*@C 159176b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 159276b2cf59SMatthew Knepley at the beginning of every step of the iteration. 159376b2cf59SMatthew Knepley 159476b2cf59SMatthew Knepley Collective on SNES 159576b2cf59SMatthew Knepley 159676b2cf59SMatthew Knepley Input Parameters: 159776b2cf59SMatthew Knepley . snes - The nonlinear solver context 159876b2cf59SMatthew Knepley . func - The function 159976b2cf59SMatthew Knepley 160076b2cf59SMatthew Knepley Calling sequence of func: 1601b5d30489SBarry Smith . func (SNES snes, PetscInt step); 160276b2cf59SMatthew Knepley 160376b2cf59SMatthew Knepley . step - The current step of the iteration 160476b2cf59SMatthew Knepley 160576b2cf59SMatthew Knepley Level: intermediate 160676b2cf59SMatthew Knepley 160776b2cf59SMatthew Knepley .keywords: SNES, update 1608b5d30489SBarry Smith 160976b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 161076b2cf59SMatthew Knepley @*/ 161177431f27SBarry Smith PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 161276b2cf59SMatthew Knepley { 161376b2cf59SMatthew Knepley PetscFunctionBegin; 16144482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 161576b2cf59SMatthew Knepley snes->update = func; 161676b2cf59SMatthew Knepley PetscFunctionReturn(0); 161776b2cf59SMatthew Knepley } 161876b2cf59SMatthew Knepley 1619e74ef692SMatthew Knepley #undef __FUNCT__ 1620e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 162176b2cf59SMatthew Knepley /*@ 162276b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 162376b2cf59SMatthew Knepley 162476b2cf59SMatthew Knepley Not collective 162576b2cf59SMatthew Knepley 162676b2cf59SMatthew Knepley Input Parameters: 162776b2cf59SMatthew Knepley . snes - The nonlinear solver context 162876b2cf59SMatthew Knepley . step - The current step of the iteration 162976b2cf59SMatthew Knepley 1630205452f4SMatthew Knepley Level: intermediate 1631205452f4SMatthew Knepley 163276b2cf59SMatthew Knepley .keywords: SNES, update 163376b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 163476b2cf59SMatthew Knepley @*/ 163577431f27SBarry Smith PetscErrorCode SNESDefaultUpdate(SNES snes, PetscInt step) 163676b2cf59SMatthew Knepley { 163776b2cf59SMatthew Knepley PetscFunctionBegin; 163876b2cf59SMatthew Knepley PetscFunctionReturn(0); 163976b2cf59SMatthew Knepley } 164076b2cf59SMatthew Knepley 16414a2ae208SSatish Balay #undef __FUNCT__ 16424a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 16439b94acceSBarry Smith /* 16449b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 16459b94acceSBarry Smith positive parameter delta. 16469b94acceSBarry Smith 16479b94acceSBarry Smith Input Parameters: 1648c7afd0dbSLois Curfman McInnes + snes - the SNES context 16499b94acceSBarry Smith . y - approximate solution of linear system 16509b94acceSBarry Smith . fnorm - 2-norm of current function 1651c7afd0dbSLois Curfman McInnes - delta - trust region size 16529b94acceSBarry Smith 16539b94acceSBarry Smith Output Parameters: 1654c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 16559b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 16569b94acceSBarry Smith region, and exceeds zero otherwise. 1657c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 16589b94acceSBarry Smith 16599b94acceSBarry Smith Note: 16604b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 16619b94acceSBarry Smith is set to be the maximum allowable step size. 16629b94acceSBarry Smith 16639b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16649b94acceSBarry Smith */ 1665dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 16669b94acceSBarry Smith { 1667064f8208SBarry Smith PetscReal nrm; 1668ea709b57SSatish Balay PetscScalar cnorm; 1669dfbe8321SBarry Smith PetscErrorCode ierr; 16703a40ed3dSBarry Smith 16713a40ed3dSBarry Smith PetscFunctionBegin; 16724482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16734482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1674c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,2); 1675184914b5SBarry Smith 1676064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1677064f8208SBarry Smith if (nrm > *delta) { 1678064f8208SBarry Smith nrm = *delta/nrm; 1679064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1680064f8208SBarry Smith cnorm = nrm; 1681329f5518SBarry Smith ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 16829b94acceSBarry Smith *ynorm = *delta; 16839b94acceSBarry Smith } else { 16849b94acceSBarry Smith *gpnorm = 0.0; 1685064f8208SBarry Smith *ynorm = nrm; 16869b94acceSBarry Smith } 16873a40ed3dSBarry Smith PetscFunctionReturn(0); 16889b94acceSBarry Smith } 16899b94acceSBarry Smith 16905968eb51SBarry Smith static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero", 16915968eb51SBarry Smith "not currently used", 16925968eb51SBarry Smith "line search failed", 16935968eb51SBarry Smith "reach maximum number of iterations", 16945968eb51SBarry Smith "function norm became NaN (not a number)", 16955968eb51SBarry Smith "not currently used", 16965968eb51SBarry Smith "number of function computations exceeded", 16975968eb51SBarry Smith "not currently used", 16985968eb51SBarry Smith "still iterating", 16995968eb51SBarry Smith "not currently used", 17005968eb51SBarry Smith "absolute size of function norm", 17015968eb51SBarry Smith "relative decrease in function norm", 17025968eb51SBarry Smith "step size is small", 17035968eb51SBarry Smith "not currently used", 17045968eb51SBarry Smith "not currently used", 17055968eb51SBarry Smith "small trust region"}; 17065968eb51SBarry Smith 17074a2ae208SSatish Balay #undef __FUNCT__ 17084a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 17099b94acceSBarry Smith /*@ 17109b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 171163a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 17129b94acceSBarry Smith 1713c7afd0dbSLois Curfman McInnes Collective on SNES 1714c7afd0dbSLois Curfman McInnes 1715b2002411SLois Curfman McInnes Input Parameters: 1716c7afd0dbSLois Curfman McInnes + snes - the SNES context 1717c7afd0dbSLois Curfman McInnes - x - the solution vector 17189b94acceSBarry Smith 1719b2002411SLois Curfman McInnes Notes: 17208ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 17218ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 17228ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 17238ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 17248ddd3da0SLois Curfman McInnes 172536851e7fSLois Curfman McInnes Level: beginner 172636851e7fSLois Curfman McInnes 17279b94acceSBarry Smith .keywords: SNES, nonlinear, solve 17289b94acceSBarry Smith 17293ab0aad5SBarry Smith .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs() 17309b94acceSBarry Smith @*/ 1731dfbe8321SBarry Smith PetscErrorCode SNESSolve(SNES snes,Vec x) 17329b94acceSBarry Smith { 1733dfbe8321SBarry Smith PetscErrorCode ierr; 1734f1af5d2fSBarry Smith PetscTruth flg; 1735052efed2SBarry Smith 17363a40ed3dSBarry Smith PetscFunctionBegin; 17374482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17384482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1739c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 17401302d50aSBarry Smith if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 174174637425SBarry Smith 174282bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1743c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 1744abc0a331SBarry Smith if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1745d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 174650ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1747c9780b6fSBarry Smith ierr = (*(snes)->solve)(snes);CHKERRQ(ierr); 1748d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1749b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 17508b179ff8SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1751da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1752da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 17535968eb51SBarry Smith if (snes->printreason) { 17545968eb51SBarry Smith if (snes->reason > 0) { 17555968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 17565968eb51SBarry Smith } else { 17575968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 17585968eb51SBarry Smith } 17595968eb51SBarry Smith } 17605968eb51SBarry Smith 17613a40ed3dSBarry Smith PetscFunctionReturn(0); 17629b94acceSBarry Smith } 17639b94acceSBarry Smith 17649b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 17659b94acceSBarry Smith 17664a2ae208SSatish Balay #undef __FUNCT__ 17674a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 176882bf6240SBarry Smith /*@C 17694b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 17709b94acceSBarry Smith 1771fee21e36SBarry Smith Collective on SNES 1772fee21e36SBarry Smith 1773c7afd0dbSLois Curfman McInnes Input Parameters: 1774c7afd0dbSLois Curfman McInnes + snes - the SNES context 1775454a90a3SBarry Smith - type - a known method 1776c7afd0dbSLois Curfman McInnes 1777c7afd0dbSLois Curfman McInnes Options Database Key: 1778454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1779c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1780ae12b187SLois Curfman McInnes 17819b94acceSBarry Smith Notes: 1782e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 17834b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1784c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17854b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1786c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17879b94acceSBarry Smith 1788ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1789ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1790ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1791ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1792ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1793ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1794ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1795ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1796ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1797b0a32e0cSBarry Smith appropriate method. 179836851e7fSLois Curfman McInnes 179936851e7fSLois Curfman McInnes Level: intermediate 1800a703fe33SLois Curfman McInnes 1801454a90a3SBarry Smith .keywords: SNES, set, type 1802435da068SBarry Smith 1803435da068SBarry Smith .seealso: SNESType, SNESCreate() 1804435da068SBarry Smith 18059b94acceSBarry Smith @*/ 1806dfbe8321SBarry Smith PetscErrorCode SNESSetType(SNES snes,const SNESType type) 18079b94acceSBarry Smith { 1808dfbe8321SBarry Smith PetscErrorCode ierr,(*r)(SNES); 18096831982aSBarry Smith PetscTruth match; 18103a40ed3dSBarry Smith 18113a40ed3dSBarry Smith PetscFunctionBegin; 18124482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18134482741eSBarry Smith PetscValidCharPointer(type,2); 181482bf6240SBarry Smith 18156831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 18160f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 181782bf6240SBarry Smith 181882bf6240SBarry Smith if (snes->setupcalled) { 1819e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 182082bf6240SBarry Smith snes->data = 0; 182182bf6240SBarry Smith } 18227d1a2b2bSBarry Smith 18239b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 182482bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1825b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1826958c9bccSBarry Smith if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1827606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 182882bf6240SBarry Smith snes->data = 0; 18293a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1830454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 18313a40ed3dSBarry Smith PetscFunctionReturn(0); 18329b94acceSBarry Smith } 18339b94acceSBarry Smith 1834a847f771SSatish Balay 18359b94acceSBarry Smith /* --------------------------------------------------------------------- */ 18364a2ae208SSatish Balay #undef __FUNCT__ 18374a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 18389b94acceSBarry Smith /*@C 18399b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1840f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 18419b94acceSBarry Smith 1842fee21e36SBarry Smith Not Collective 1843fee21e36SBarry Smith 184436851e7fSLois Curfman McInnes Level: advanced 184536851e7fSLois Curfman McInnes 18469b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 18479b94acceSBarry Smith 18489b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 18499b94acceSBarry Smith @*/ 1850dfbe8321SBarry Smith PetscErrorCode SNESRegisterDestroy(void) 18519b94acceSBarry Smith { 1852dfbe8321SBarry Smith PetscErrorCode ierr; 185382bf6240SBarry Smith 18543a40ed3dSBarry Smith PetscFunctionBegin; 185582bf6240SBarry Smith if (SNESList) { 1856b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 185782bf6240SBarry Smith SNESList = 0; 18589b94acceSBarry Smith } 18594c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 18603a40ed3dSBarry Smith PetscFunctionReturn(0); 18619b94acceSBarry Smith } 18629b94acceSBarry Smith 18634a2ae208SSatish Balay #undef __FUNCT__ 18644a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 18659b94acceSBarry Smith /*@C 18669a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 18679b94acceSBarry Smith 1868c7afd0dbSLois Curfman McInnes Not Collective 1869c7afd0dbSLois Curfman McInnes 18709b94acceSBarry Smith Input Parameter: 18714b0e389bSBarry Smith . snes - nonlinear solver context 18729b94acceSBarry Smith 18739b94acceSBarry Smith Output Parameter: 18743a7fca6bSBarry Smith . type - SNES method (a character string) 18759b94acceSBarry Smith 187636851e7fSLois Curfman McInnes Level: intermediate 187736851e7fSLois Curfman McInnes 1878454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 18799b94acceSBarry Smith @*/ 1880dfbe8321SBarry Smith PetscErrorCode SNESGetType(SNES snes,SNESType *type) 18819b94acceSBarry Smith { 18823a40ed3dSBarry Smith PetscFunctionBegin; 18834482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18844482741eSBarry Smith PetscValidPointer(type,2); 1885454a90a3SBarry Smith *type = snes->type_name; 18863a40ed3dSBarry Smith PetscFunctionReturn(0); 18879b94acceSBarry Smith } 18889b94acceSBarry Smith 18894a2ae208SSatish Balay #undef __FUNCT__ 18904a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 18919b94acceSBarry Smith /*@C 18929b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 18939b94acceSBarry Smith stored. 18949b94acceSBarry Smith 1895c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1896c7afd0dbSLois Curfman McInnes 18979b94acceSBarry Smith Input Parameter: 18989b94acceSBarry Smith . snes - the SNES context 18999b94acceSBarry Smith 19009b94acceSBarry Smith Output Parameter: 19019b94acceSBarry Smith . x - the solution 19029b94acceSBarry Smith 190336851e7fSLois Curfman McInnes Level: advanced 190436851e7fSLois Curfman McInnes 19059b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 19069b94acceSBarry Smith 19074b27c08aSLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetSolutionUpdate() 19089b94acceSBarry Smith @*/ 1909dfbe8321SBarry Smith PetscErrorCode SNESGetSolution(SNES snes,Vec *x) 19109b94acceSBarry Smith { 19113a40ed3dSBarry Smith PetscFunctionBegin; 19124482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19134482741eSBarry Smith PetscValidPointer(x,2); 19149b94acceSBarry Smith *x = snes->vec_sol_always; 19153a40ed3dSBarry Smith PetscFunctionReturn(0); 19169b94acceSBarry Smith } 19179b94acceSBarry Smith 19184a2ae208SSatish Balay #undef __FUNCT__ 19194a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 19209b94acceSBarry Smith /*@C 19219b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 19229b94acceSBarry Smith stored. 19239b94acceSBarry Smith 1924c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1925c7afd0dbSLois Curfman McInnes 19269b94acceSBarry Smith Input Parameter: 19279b94acceSBarry Smith . snes - the SNES context 19289b94acceSBarry Smith 19299b94acceSBarry Smith Output Parameter: 19309b94acceSBarry Smith . x - the solution update 19319b94acceSBarry Smith 193236851e7fSLois Curfman McInnes Level: advanced 193336851e7fSLois Curfman McInnes 19349b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 19359b94acceSBarry Smith 19369b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 19379b94acceSBarry Smith @*/ 1938dfbe8321SBarry Smith PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x) 19399b94acceSBarry Smith { 19403a40ed3dSBarry Smith PetscFunctionBegin; 19414482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19424482741eSBarry Smith PetscValidPointer(x,2); 19439b94acceSBarry Smith *x = snes->vec_sol_update_always; 19443a40ed3dSBarry Smith PetscFunctionReturn(0); 19459b94acceSBarry Smith } 19469b94acceSBarry Smith 19474a2ae208SSatish Balay #undef __FUNCT__ 19484a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 19499b94acceSBarry Smith /*@C 19503638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 19519b94acceSBarry Smith 1952c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1953c7afd0dbSLois Curfman McInnes 19549b94acceSBarry Smith Input Parameter: 19559b94acceSBarry Smith . snes - the SNES context 19569b94acceSBarry Smith 19579b94acceSBarry Smith Output Parameter: 19587bf4e008SBarry Smith + r - the function (or PETSC_NULL) 195900036973SBarry Smith . ctx - the function context (or PETSC_NULL) 196000036973SBarry Smith - func - the function (or PETSC_NULL) 19619b94acceSBarry Smith 196236851e7fSLois Curfman McInnes Level: advanced 196336851e7fSLois Curfman McInnes 1964a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 19659b94acceSBarry Smith 19664b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 19679b94acceSBarry Smith @*/ 19686849ba73SBarry Smith PetscErrorCode SNESGetFunction(SNES snes,Vec *r,void **ctx,PetscErrorCode (**func)(SNES,Vec,Vec,void*)) 19699b94acceSBarry Smith { 19703a40ed3dSBarry Smith PetscFunctionBegin; 19714482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19727bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 19737bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 197400036973SBarry Smith if (func) *func = snes->computefunction; 19753a40ed3dSBarry Smith PetscFunctionReturn(0); 19769b94acceSBarry Smith } 19779b94acceSBarry Smith 19784a2ae208SSatish Balay #undef __FUNCT__ 19794a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 19803c7409f5SSatish Balay /*@C 19813c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1982d850072dSLois Curfman McInnes SNES options in the database. 19833c7409f5SSatish Balay 1984fee21e36SBarry Smith Collective on SNES 1985fee21e36SBarry Smith 1986c7afd0dbSLois Curfman McInnes Input Parameter: 1987c7afd0dbSLois Curfman McInnes + snes - the SNES context 1988c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1989c7afd0dbSLois Curfman McInnes 1990d850072dSLois Curfman McInnes Notes: 1991a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1992c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1993d850072dSLois Curfman McInnes 199436851e7fSLois Curfman McInnes Level: advanced 199536851e7fSLois Curfman McInnes 19963c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1997a86d99e1SLois Curfman McInnes 1998a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19993c7409f5SSatish Balay @*/ 2000dfbe8321SBarry Smith PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[]) 20013c7409f5SSatish Balay { 2002dfbe8321SBarry Smith PetscErrorCode ierr; 20033c7409f5SSatish Balay 20043a40ed3dSBarry Smith PetscFunctionBegin; 20054482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2006639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 200794b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20083a40ed3dSBarry Smith PetscFunctionReturn(0); 20093c7409f5SSatish Balay } 20103c7409f5SSatish Balay 20114a2ae208SSatish Balay #undef __FUNCT__ 20124a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 20133c7409f5SSatish Balay /*@C 2014f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2015d850072dSLois Curfman McInnes SNES options in the database. 20163c7409f5SSatish Balay 2017fee21e36SBarry Smith Collective on SNES 2018fee21e36SBarry Smith 2019c7afd0dbSLois Curfman McInnes Input Parameters: 2020c7afd0dbSLois Curfman McInnes + snes - the SNES context 2021c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2022c7afd0dbSLois Curfman McInnes 2023d850072dSLois Curfman McInnes Notes: 2024a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2025c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2026d850072dSLois Curfman McInnes 202736851e7fSLois Curfman McInnes Level: advanced 202836851e7fSLois Curfman McInnes 20293c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2030a86d99e1SLois Curfman McInnes 2031a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20323c7409f5SSatish Balay @*/ 2033dfbe8321SBarry Smith PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 20343c7409f5SSatish Balay { 2035dfbe8321SBarry Smith PetscErrorCode ierr; 20363c7409f5SSatish Balay 20373a40ed3dSBarry Smith PetscFunctionBegin; 20384482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2039639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 204094b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20413a40ed3dSBarry Smith PetscFunctionReturn(0); 20423c7409f5SSatish Balay } 20433c7409f5SSatish Balay 20444a2ae208SSatish Balay #undef __FUNCT__ 20454a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 20469ab63eb5SSatish Balay /*@C 20473c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20483c7409f5SSatish Balay SNES options in the database. 20493c7409f5SSatish Balay 2050c7afd0dbSLois Curfman McInnes Not Collective 2051c7afd0dbSLois Curfman McInnes 20523c7409f5SSatish Balay Input Parameter: 20533c7409f5SSatish Balay . snes - the SNES context 20543c7409f5SSatish Balay 20553c7409f5SSatish Balay Output Parameter: 20563c7409f5SSatish Balay . prefix - pointer to the prefix string used 20573c7409f5SSatish Balay 20589ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 20599ab63eb5SSatish Balay sufficient length to hold the prefix. 20609ab63eb5SSatish Balay 206136851e7fSLois Curfman McInnes Level: advanced 206236851e7fSLois Curfman McInnes 20633c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2064a86d99e1SLois Curfman McInnes 2065a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20663c7409f5SSatish Balay @*/ 2067dfbe8321SBarry Smith PetscErrorCode SNESGetOptionsPrefix(SNES snes,char *prefix[]) 20683c7409f5SSatish Balay { 2069dfbe8321SBarry Smith PetscErrorCode ierr; 20703c7409f5SSatish Balay 20713a40ed3dSBarry Smith PetscFunctionBegin; 20724482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2073639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 20743a40ed3dSBarry Smith PetscFunctionReturn(0); 20753c7409f5SSatish Balay } 20763c7409f5SSatish Balay 2077b2002411SLois Curfman McInnes 20784a2ae208SSatish Balay #undef __FUNCT__ 20794a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 20803cea93caSBarry Smith /*@C 20813cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 20823cea93caSBarry Smith 20837f6c08e0SMatthew Knepley Level: advanced 20843cea93caSBarry Smith @*/ 20856849ba73SBarry Smith PetscErrorCode SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2086b2002411SLois Curfman McInnes { 2087e2d1d2b7SBarry Smith char fullname[PETSC_MAX_PATH_LEN]; 2088dfbe8321SBarry Smith PetscErrorCode ierr; 2089b2002411SLois Curfman McInnes 2090b2002411SLois Curfman McInnes PetscFunctionBegin; 2091b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2092c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2093b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2094b2002411SLois Curfman McInnes } 2095da9b6338SBarry Smith 2096da9b6338SBarry Smith #undef __FUNCT__ 2097da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 2098dfbe8321SBarry Smith PetscErrorCode SNESTestLocalMin(SNES snes) 2099da9b6338SBarry Smith { 2100dfbe8321SBarry Smith PetscErrorCode ierr; 210177431f27SBarry Smith PetscInt N,i,j; 2102da9b6338SBarry Smith Vec u,uh,fh; 2103da9b6338SBarry Smith PetscScalar value; 2104da9b6338SBarry Smith PetscReal norm; 2105da9b6338SBarry Smith 2106da9b6338SBarry Smith PetscFunctionBegin; 2107da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2108da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2109da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2110da9b6338SBarry Smith 2111da9b6338SBarry Smith /* currently only works for sequential */ 2112da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2113da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2114da9b6338SBarry Smith for (i=0; i<N; i++) { 2115da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 211677431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2117da9b6338SBarry Smith for (j=-10; j<11; j++) { 2118ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2119da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 21203ab0aad5SBarry Smith ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2121da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 212277431f27SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2123da9b6338SBarry Smith value = -value; 2124da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2125da9b6338SBarry Smith } 2126da9b6338SBarry Smith } 2127da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2128da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2129da9b6338SBarry Smith PetscFunctionReturn(0); 2130da9b6338SBarry Smith } 2131