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 } 75b0a32e0cSBarry 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", 7777d8c4bbSBarry Smith snes->rtol,snes->atol,snes->xtol);CHKERRQ(ierr); 78b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr); 79b0a32e0cSBarry 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) { 83b0a32e0cSBarry 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 10476b2cf59SMatthew Knepley static int 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) { 12576b2cf59SMatthew Knepley SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %d allowed", MAXSETFROMOPTIONS); 12676b2cf59SMatthew Knepley } 12776b2cf59SMatthew Knepley 12876b2cf59SMatthew Knepley othersetfromoptions[numberofsetfromoptions++] = snescheck; 12976b2cf59SMatthew Knepley PetscFunctionReturn(0); 13076b2cf59SMatthew Knepley } 13176b2cf59SMatthew Knepley 1324a2ae208SSatish Balay #undef __FUNCT__ 1334a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions" 1349b94acceSBarry Smith /*@ 13594b7f48cSBarry Smith SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 1369b94acceSBarry Smith 137c7afd0dbSLois Curfman McInnes Collective on SNES 138c7afd0dbSLois Curfman McInnes 1399b94acceSBarry Smith Input Parameter: 1409b94acceSBarry Smith . snes - the SNES context 1419b94acceSBarry Smith 14236851e7fSLois Curfman McInnes Options Database Keys: 1436831982aSBarry Smith + -snes_type <type> - ls, tr, umls, umtr, test 14482738288SBarry Smith . -snes_stol - convergence tolerance in terms of the norm 14582738288SBarry Smith of the change in the solution between steps 146b39c3a46SLois Curfman McInnes . -snes_atol <atol> - absolute tolerance of residual norm 147b39c3a46SLois Curfman McInnes . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 148b39c3a46SLois Curfman McInnes . -snes_max_it <max_it> - maximum number of iterations 149b39c3a46SLois Curfman McInnes . -snes_max_funcs <max_funcs> - maximum number of function evaluations 15050ffb88aSMatthew Knepley . -snes_max_fail <max_fail> - maximum number of failures 151b39c3a46SLois Curfman McInnes . -snes_trtol <trtol> - trust region tolerance 15282738288SBarry Smith . -snes_no_convergence_test - skip convergence test in nonlinear or minimization 15382738288SBarry Smith solver; hence iterations will continue until max_it 1541fbbfb26SLois Curfman McInnes or some other criterion is reached. Saves expense 15582738288SBarry Smith of convergence test 15682738288SBarry Smith . -snes_monitor - prints residual norm at each iteration 157d132466eSBarry Smith . -snes_vecmonitor - plots solution at each iteration 158d132466eSBarry Smith . -snes_vecmonitor_update - plots update to solution at each iteration 15982738288SBarry Smith . -snes_xmonitor - plots residual norm at each iteration 160e24b481bSBarry Smith . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 1615968eb51SBarry Smith . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 1625968eb51SBarry Smith - -snes_print_converged_reason - print the reason for convergence/divergence after each solve 16382738288SBarry Smith 16482738288SBarry Smith Options Database for Eisenstat-Walker method: 1654b27c08aSLois Curfman McInnes + -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence 1664b27c08aSLois Curfman McInnes . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 16736851e7fSLois Curfman McInnes . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 16836851e7fSLois Curfman McInnes . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 16936851e7fSLois Curfman McInnes . -snes_ksp_ew_gamma <gamma> - Sets gamma 17036851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha <alpha> - Sets alpha 17136851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 17236851e7fSLois Curfman McInnes - -snes_ksp_ew_threshold <threshold> - Sets threshold 17382738288SBarry Smith 17411ca99fdSLois Curfman McInnes Notes: 17511ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 17611ca99fdSLois Curfman McInnes the users manual. 17783e2fdc7SBarry Smith 17836851e7fSLois Curfman McInnes Level: beginner 17936851e7fSLois Curfman McInnes 1809b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1819b94acceSBarry Smith 18269ed319cSSatish Balay .seealso: SNESSetOptionsPrefix() 1839b94acceSBarry Smith @*/ 184dfbe8321SBarry Smith PetscErrorCode SNESSetFromOptions(SNES snes) 1859b94acceSBarry Smith { 18694b7f48cSBarry Smith KSP ksp; 187186905e3SBarry Smith SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 188f1af5d2fSBarry Smith PetscTruth flg; 189dfbe8321SBarry Smith PetscErrorCode ierr; 190dfbe8321SBarry Smith int i; 1912fc52814SBarry Smith const char *deft; 1922fc52814SBarry Smith char type[256]; 1939b94acceSBarry Smith 1943a40ed3dSBarry Smith PetscFunctionBegin; 1954482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 196ca161407SBarry Smith 197b0a32e0cSBarry Smith ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 198186905e3SBarry Smith if (snes->type_name) { 199186905e3SBarry Smith deft = snes->type_name; 200186905e3SBarry Smith } else { 2014b27c08aSLois Curfman McInnes deft = SNESLS; 202d64ed03dSBarry Smith } 2034bbc92c1SBarry Smith 204186905e3SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 205b0a32e0cSBarry Smith ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 206d64ed03dSBarry Smith if (flg) { 207186905e3SBarry Smith ierr = SNESSetType(snes,type);CHKERRQ(ierr); 208186905e3SBarry Smith } else if (!snes->type_name) { 209186905e3SBarry Smith ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 210d64ed03dSBarry Smith } 211909c8a9fSBarry Smith ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 21293c39befSBarry Smith 21387828ca2SBarry Smith ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 21487828ca2SBarry Smith ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr); 215186905e3SBarry Smith 21687828ca2SBarry Smith ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 217b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 218b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 21950ffb88aSMatthew Knepley ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 2205968eb51SBarry Smith ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr); 2215968eb51SBarry Smith if (flg) { 2225968eb51SBarry Smith snes->printreason = PETSC_TRUE; 2235968eb51SBarry Smith } 224186905e3SBarry Smith 225b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr); 226186905e3SBarry Smith 227b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 22887828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 22987828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 23087828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 23187828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 23287828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 23387828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 234186905e3SBarry Smith 235b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr); 23693c39befSBarry Smith if (flg) {snes->converged = 0;} 237b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr); 2385cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 239b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr); 240b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);} 2413a7fca6bSBarry Smith ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr); 2423a7fca6bSBarry Smith if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);} 243b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr); 244b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);} 245b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr); 246b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 247b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr); 2487c922b88SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);} 2495ed2d596SBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr); 2505ed2d596SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);} 251b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr); 252186905e3SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 253e24b481bSBarry Smith 254b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 2554b27c08aSLois Curfman McInnes if (flg) { 256186905e3SBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 257b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 2589b94acceSBarry Smith } 259639f9d9dSBarry Smith 26076b2cf59SMatthew Knepley for(i = 0; i < numberofsetfromoptions; i++) { 26176b2cf59SMatthew Knepley ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 26276b2cf59SMatthew Knepley } 26376b2cf59SMatthew Knepley 264186905e3SBarry Smith if (snes->setfromoptions) { 265186905e3SBarry Smith ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 266639f9d9dSBarry Smith } 267639f9d9dSBarry Smith 268b0a32e0cSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 2694bbc92c1SBarry Smith 27094b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 27194b7f48cSBarry Smith ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 27293993e2dSLois Curfman McInnes 2733a40ed3dSBarry Smith PetscFunctionReturn(0); 2749b94acceSBarry Smith } 2759b94acceSBarry Smith 276a847f771SSatish Balay 2774a2ae208SSatish Balay #undef __FUNCT__ 2784a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext" 2799b94acceSBarry Smith /*@ 2809b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 2819b94acceSBarry Smith the nonlinear solvers. 2829b94acceSBarry Smith 283fee21e36SBarry Smith Collective on SNES 284fee21e36SBarry Smith 285c7afd0dbSLois Curfman McInnes Input Parameters: 286c7afd0dbSLois Curfman McInnes + snes - the SNES context 287c7afd0dbSLois Curfman McInnes - usrP - optional user context 288c7afd0dbSLois Curfman McInnes 28936851e7fSLois Curfman McInnes Level: intermediate 29036851e7fSLois Curfman McInnes 2919b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 2929b94acceSBarry Smith 2939b94acceSBarry Smith .seealso: SNESGetApplicationContext() 2949b94acceSBarry Smith @*/ 295dfbe8321SBarry Smith PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP) 2969b94acceSBarry Smith { 2973a40ed3dSBarry Smith PetscFunctionBegin; 2984482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2999b94acceSBarry Smith snes->user = usrP; 3003a40ed3dSBarry Smith PetscFunctionReturn(0); 3019b94acceSBarry Smith } 30274679c65SBarry Smith 3034a2ae208SSatish Balay #undef __FUNCT__ 3044a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext" 3059b94acceSBarry Smith /*@C 3069b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3079b94acceSBarry Smith nonlinear solvers. 3089b94acceSBarry Smith 309c7afd0dbSLois Curfman McInnes Not Collective 310c7afd0dbSLois Curfman McInnes 3119b94acceSBarry Smith Input Parameter: 3129b94acceSBarry Smith . snes - SNES context 3139b94acceSBarry Smith 3149b94acceSBarry Smith Output Parameter: 3159b94acceSBarry Smith . usrP - user context 3169b94acceSBarry Smith 31736851e7fSLois Curfman McInnes Level: intermediate 31836851e7fSLois Curfman McInnes 3199b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3209b94acceSBarry Smith 3219b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3229b94acceSBarry Smith @*/ 323dfbe8321SBarry Smith PetscErrorCode SNESGetApplicationContext(SNES snes,void **usrP) 3249b94acceSBarry Smith { 3253a40ed3dSBarry Smith PetscFunctionBegin; 3264482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 3279b94acceSBarry Smith *usrP = snes->user; 3283a40ed3dSBarry Smith PetscFunctionReturn(0); 3299b94acceSBarry Smith } 33074679c65SBarry Smith 3314a2ae208SSatish Balay #undef __FUNCT__ 3324a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber" 3339b94acceSBarry Smith /*@ 334c8228a4eSBarry Smith SNESGetIterationNumber - Gets the number of nonlinear iterations completed 335c8228a4eSBarry Smith at this time. 3369b94acceSBarry Smith 337c7afd0dbSLois Curfman McInnes Not Collective 338c7afd0dbSLois Curfman McInnes 3399b94acceSBarry Smith Input Parameter: 3409b94acceSBarry Smith . snes - SNES context 3419b94acceSBarry Smith 3429b94acceSBarry Smith Output Parameter: 3439b94acceSBarry Smith . iter - iteration number 3449b94acceSBarry Smith 345c8228a4eSBarry Smith Notes: 346c8228a4eSBarry Smith For example, during the computation of iteration 2 this would return 1. 347c8228a4eSBarry Smith 348c8228a4eSBarry Smith This is useful for using lagged Jacobians (where one does not recompute the 34908405cd6SLois Curfman McInnes Jacobian at each SNES iteration). For example, the code 35008405cd6SLois Curfman McInnes .vb 35108405cd6SLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&it); 35208405cd6SLois Curfman McInnes if (!(it % 2)) { 35308405cd6SLois Curfman McInnes [compute Jacobian here] 35408405cd6SLois Curfman McInnes } 35508405cd6SLois Curfman McInnes .ve 356c8228a4eSBarry Smith can be used in your ComputeJacobian() function to cause the Jacobian to be 35708405cd6SLois Curfman McInnes recomputed every second SNES iteration. 358c8228a4eSBarry Smith 35936851e7fSLois Curfman McInnes Level: intermediate 36036851e7fSLois Curfman McInnes 3619b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3629b94acceSBarry Smith @*/ 363dfbe8321SBarry Smith PetscErrorCode SNESGetIterationNumber(SNES snes,int* iter) 3649b94acceSBarry Smith { 3653a40ed3dSBarry Smith PetscFunctionBegin; 3664482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 3674482741eSBarry Smith PetscValidIntPointer(iter,2); 3689b94acceSBarry Smith *iter = snes->iter; 3693a40ed3dSBarry Smith PetscFunctionReturn(0); 3709b94acceSBarry Smith } 37174679c65SBarry Smith 3724a2ae208SSatish Balay #undef __FUNCT__ 3734a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm" 3749b94acceSBarry Smith /*@ 3759b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3769b94acceSBarry Smith with SNESSSetFunction(). 3779b94acceSBarry Smith 378c7afd0dbSLois Curfman McInnes Collective on SNES 379c7afd0dbSLois Curfman McInnes 3809b94acceSBarry Smith Input Parameter: 3819b94acceSBarry Smith . snes - SNES context 3829b94acceSBarry Smith 3839b94acceSBarry Smith Output Parameter: 3849b94acceSBarry Smith . fnorm - 2-norm of function 3859b94acceSBarry Smith 38636851e7fSLois Curfman McInnes Level: intermediate 38736851e7fSLois Curfman McInnes 3889b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 389a86d99e1SLois Curfman McInnes 39008405cd6SLois Curfman McInnes .seealso: SNESGetFunction() 3919b94acceSBarry Smith @*/ 392dfbe8321SBarry Smith PetscErrorCode SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm) 3939b94acceSBarry Smith { 3943a40ed3dSBarry Smith PetscFunctionBegin; 3954482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 3964482741eSBarry Smith PetscValidScalarPointer(fnorm,2); 3979b94acceSBarry Smith *fnorm = snes->norm; 3983a40ed3dSBarry Smith PetscFunctionReturn(0); 3999b94acceSBarry Smith } 40074679c65SBarry Smith 4014a2ae208SSatish Balay #undef __FUNCT__ 4024a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps" 4039b94acceSBarry Smith /*@ 4049b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4059b94acceSBarry Smith attempted by the nonlinear solver. 4069b94acceSBarry Smith 407c7afd0dbSLois Curfman McInnes Not Collective 408c7afd0dbSLois Curfman McInnes 4099b94acceSBarry Smith Input Parameter: 4109b94acceSBarry Smith . snes - SNES context 4119b94acceSBarry Smith 4129b94acceSBarry Smith Output Parameter: 4139b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4149b94acceSBarry Smith 415c96a6f78SLois Curfman McInnes Notes: 416c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 417c96a6f78SLois Curfman McInnes 41836851e7fSLois Curfman McInnes Level: intermediate 41936851e7fSLois Curfman McInnes 4209b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4219b94acceSBarry Smith @*/ 422dfbe8321SBarry Smith PetscErrorCode SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4239b94acceSBarry Smith { 4243a40ed3dSBarry Smith PetscFunctionBegin; 4254482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 4264482741eSBarry Smith PetscValidIntPointer(nfails,2); 42750ffb88aSMatthew Knepley *nfails = snes->numFailures; 42850ffb88aSMatthew Knepley PetscFunctionReturn(0); 42950ffb88aSMatthew Knepley } 43050ffb88aSMatthew Knepley 43150ffb88aSMatthew Knepley #undef __FUNCT__ 43250ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps" 43350ffb88aSMatthew Knepley /*@ 43450ffb88aSMatthew Knepley SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps 43550ffb88aSMatthew Knepley attempted by the nonlinear solver before it gives up. 43650ffb88aSMatthew Knepley 43750ffb88aSMatthew Knepley Not Collective 43850ffb88aSMatthew Knepley 43950ffb88aSMatthew Knepley Input Parameters: 44050ffb88aSMatthew Knepley + snes - SNES context 44150ffb88aSMatthew Knepley - maxFails - maximum of unsuccessful steps 44250ffb88aSMatthew Knepley 44350ffb88aSMatthew Knepley Level: intermediate 44450ffb88aSMatthew Knepley 44550ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 44650ffb88aSMatthew Knepley @*/ 447dfbe8321SBarry Smith PetscErrorCode SNESSetMaximumUnsuccessfulSteps(SNES snes, int maxFails) 44850ffb88aSMatthew Knepley { 44950ffb88aSMatthew Knepley PetscFunctionBegin; 4504482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 45150ffb88aSMatthew Knepley snes->maxFailures = maxFails; 45250ffb88aSMatthew Knepley PetscFunctionReturn(0); 45350ffb88aSMatthew Knepley } 45450ffb88aSMatthew Knepley 45550ffb88aSMatthew Knepley #undef __FUNCT__ 45650ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps" 45750ffb88aSMatthew Knepley /*@ 45850ffb88aSMatthew Knepley SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps 45950ffb88aSMatthew Knepley attempted by the nonlinear solver before it gives up. 46050ffb88aSMatthew Knepley 46150ffb88aSMatthew Knepley Not Collective 46250ffb88aSMatthew Knepley 46350ffb88aSMatthew Knepley Input Parameter: 46450ffb88aSMatthew Knepley . snes - SNES context 46550ffb88aSMatthew Knepley 46650ffb88aSMatthew Knepley Output Parameter: 46750ffb88aSMatthew Knepley . maxFails - maximum of unsuccessful steps 46850ffb88aSMatthew Knepley 46950ffb88aSMatthew Knepley Level: intermediate 47050ffb88aSMatthew Knepley 47150ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 47250ffb88aSMatthew Knepley @*/ 473dfbe8321SBarry Smith PetscErrorCode SNESGetMaximumUnsuccessfulSteps(SNES snes, int *maxFails) 47450ffb88aSMatthew Knepley { 47550ffb88aSMatthew Knepley PetscFunctionBegin; 4764482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 4774482741eSBarry Smith PetscValidIntPointer(maxFails,2); 47850ffb88aSMatthew Knepley *maxFails = snes->maxFailures; 4793a40ed3dSBarry Smith PetscFunctionReturn(0); 4809b94acceSBarry Smith } 481a847f771SSatish Balay 4824a2ae208SSatish Balay #undef __FUNCT__ 4834a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations" 484c96a6f78SLois Curfman McInnes /*@ 485c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 486c96a6f78SLois Curfman McInnes used by the nonlinear solver. 487c96a6f78SLois Curfman McInnes 488c7afd0dbSLois Curfman McInnes Not Collective 489c7afd0dbSLois Curfman McInnes 490c96a6f78SLois Curfman McInnes Input Parameter: 491c96a6f78SLois Curfman McInnes . snes - SNES context 492c96a6f78SLois Curfman McInnes 493c96a6f78SLois Curfman McInnes Output Parameter: 494c96a6f78SLois Curfman McInnes . lits - number of linear iterations 495c96a6f78SLois Curfman McInnes 496c96a6f78SLois Curfman McInnes Notes: 497c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 498c96a6f78SLois Curfman McInnes 49936851e7fSLois Curfman McInnes Level: intermediate 50036851e7fSLois Curfman McInnes 501c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 502c96a6f78SLois Curfman McInnes @*/ 503dfbe8321SBarry Smith PetscErrorCode SNESGetNumberLinearIterations(SNES snes,int* lits) 504c96a6f78SLois Curfman McInnes { 5053a40ed3dSBarry Smith PetscFunctionBegin; 5064482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 5074482741eSBarry Smith PetscValidIntPointer(lits,2); 508c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 5093a40ed3dSBarry Smith PetscFunctionReturn(0); 510c96a6f78SLois Curfman McInnes } 511c96a6f78SLois Curfman McInnes 5124a2ae208SSatish Balay #undef __FUNCT__ 51394b7f48cSBarry Smith #define __FUNCT__ "SNESGetKSP" 5149b94acceSBarry Smith /*@C 51594b7f48cSBarry Smith SNESGetKSP - Returns the KSP context for a SNES solver. 5169b94acceSBarry Smith 51794b7f48cSBarry Smith Not Collective, but if SNES object is parallel, then KSP object is parallel 518c7afd0dbSLois Curfman McInnes 5199b94acceSBarry Smith Input Parameter: 5209b94acceSBarry Smith . snes - the SNES context 5219b94acceSBarry Smith 5229b94acceSBarry Smith Output Parameter: 52394b7f48cSBarry Smith . ksp - the KSP context 5249b94acceSBarry Smith 5259b94acceSBarry Smith Notes: 52694b7f48cSBarry Smith The user can then directly manipulate the KSP context to set various 5279b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5289b94acceSBarry Smith KSP and PC contexts as well. 5299b94acceSBarry Smith 53036851e7fSLois Curfman McInnes Level: beginner 53136851e7fSLois Curfman McInnes 53294b7f48cSBarry Smith .keywords: SNES, nonlinear, get, KSP, context 5339b94acceSBarry Smith 53494b7f48cSBarry Smith .seealso: KSPGetPC() 5359b94acceSBarry Smith @*/ 536dfbe8321SBarry Smith PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp) 5379b94acceSBarry Smith { 5383a40ed3dSBarry Smith PetscFunctionBegin; 5394482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 5404482741eSBarry Smith PetscValidPointer(ksp,2); 54194b7f48cSBarry Smith *ksp = snes->ksp; 5423a40ed3dSBarry Smith PetscFunctionReturn(0); 5439b94acceSBarry Smith } 54482bf6240SBarry Smith 5454a2ae208SSatish Balay #undef __FUNCT__ 5464a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc" 5476849ba73SBarry Smith static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 548e24b481bSBarry Smith { 549e24b481bSBarry Smith PetscFunctionBegin; 550e24b481bSBarry Smith PetscFunctionReturn(0); 551e24b481bSBarry Smith } 552e24b481bSBarry Smith 5539b94acceSBarry Smith /* -----------------------------------------------------------*/ 5544a2ae208SSatish Balay #undef __FUNCT__ 5554a2ae208SSatish Balay #define __FUNCT__ "SNESCreate" 5569b94acceSBarry Smith /*@C 5579b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5589b94acceSBarry Smith 559c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 560c7afd0dbSLois Curfman McInnes 561c7afd0dbSLois Curfman McInnes Input Parameters: 562c7afd0dbSLois Curfman McInnes + comm - MPI communicator 5639b94acceSBarry Smith 5649b94acceSBarry Smith Output Parameter: 5659b94acceSBarry Smith . outsnes - the new SNES context 5669b94acceSBarry Smith 567c7afd0dbSLois Curfman McInnes Options Database Keys: 568c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 569c7afd0dbSLois Curfman McInnes and no preconditioning matrix 570c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 571c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 572c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 573c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 574c1f60f51SBarry Smith 57536851e7fSLois Curfman McInnes Level: beginner 57636851e7fSLois Curfman McInnes 5779b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5789b94acceSBarry Smith 5794b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES 5809b94acceSBarry Smith @*/ 581dfbe8321SBarry Smith PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes) 5829b94acceSBarry Smith { 583dfbe8321SBarry Smith PetscErrorCode ierr; 5849b94acceSBarry Smith SNES snes; 5859b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 58637fcc0dbSBarry Smith 5873a40ed3dSBarry Smith PetscFunctionBegin; 5884482741eSBarry Smith PetscValidPointer(outsnes,1); 5898ba1e511SMatthew Knepley *outsnes = PETSC_NULL; 5908ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 5918ba1e511SMatthew Knepley ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 5928ba1e511SMatthew Knepley #endif 5938ba1e511SMatthew Knepley 5943f1db9ecSBarry Smith PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView); 595b0a32e0cSBarry Smith PetscLogObjectCreate(snes); 596e24b481bSBarry Smith snes->bops->publish = SNESPublish_Petsc; 5979b94acceSBarry Smith snes->max_its = 50; 5989750a799SBarry Smith snes->max_funcs = 10000; 5999b94acceSBarry Smith snes->norm = 0.0; 600b4874afaSBarry Smith snes->rtol = 1.e-8; 601b4874afaSBarry Smith snes->ttol = 0.0; 602b4874afaSBarry Smith snes->atol = 1.e-50; 6039b94acceSBarry Smith snes->xtol = 1.e-8; 6044b27c08aSLois Curfman McInnes snes->deltatol = 1.e-12; 6059b94acceSBarry Smith snes->nfuncs = 0; 60650ffb88aSMatthew Knepley snes->numFailures = 0; 60750ffb88aSMatthew Knepley snes->maxFailures = 1; 6087a00f4a9SLois Curfman McInnes snes->linear_its = 0; 609639f9d9dSBarry Smith snes->numbermonitors = 0; 6109b94acceSBarry Smith snes->data = 0; 6119b94acceSBarry Smith snes->view = 0; 61282bf6240SBarry Smith snes->setupcalled = 0; 613186905e3SBarry Smith snes->ksp_ewconv = PETSC_FALSE; 6146f24a144SLois Curfman McInnes snes->vwork = 0; 6156f24a144SLois Curfman McInnes snes->nwork = 0; 616758f92a0SBarry Smith snes->conv_hist_len = 0; 617758f92a0SBarry Smith snes->conv_hist_max = 0; 618758f92a0SBarry Smith snes->conv_hist = PETSC_NULL; 619758f92a0SBarry Smith snes->conv_hist_its = PETSC_NULL; 620758f92a0SBarry Smith snes->conv_hist_reset = PETSC_TRUE; 621184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 6229b94acceSBarry Smith 6239b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 624b0a32e0cSBarry Smith ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 625b0a32e0cSBarry Smith PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 6269b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6279b94acceSBarry Smith kctx->version = 2; 6289b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6299b94acceSBarry Smith this was too large for some test cases */ 6309b94acceSBarry Smith kctx->rtol_last = 0; 6319b94acceSBarry Smith kctx->rtol_max = .9; 6329b94acceSBarry Smith kctx->gamma = 1.0; 6339b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6349b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6359b94acceSBarry Smith kctx->threshold = .1; 6369b94acceSBarry Smith kctx->lresid_last = 0; 6379b94acceSBarry Smith kctx->norm_last = 0; 6389b94acceSBarry Smith 63994b7f48cSBarry Smith ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr); 64094b7f48cSBarry Smith PetscLogObjectParent(snes,snes->ksp) 6415334005bSBarry Smith 6429b94acceSBarry Smith *outsnes = snes; 64300036973SBarry Smith ierr = PetscPublishAll(snes);CHKERRQ(ierr); 6443a40ed3dSBarry Smith PetscFunctionReturn(0); 6459b94acceSBarry Smith } 6469b94acceSBarry Smith 6474a2ae208SSatish Balay #undef __FUNCT__ 6484a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction" 6499b94acceSBarry Smith /*@C 6509b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6519b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6529b94acceSBarry Smith equations. 6539b94acceSBarry Smith 654fee21e36SBarry Smith Collective on SNES 655fee21e36SBarry Smith 656c7afd0dbSLois Curfman McInnes Input Parameters: 657c7afd0dbSLois Curfman McInnes + snes - the SNES context 658c7afd0dbSLois Curfman McInnes . func - function evaluation routine 659c7afd0dbSLois Curfman McInnes . r - vector to store function value 660c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 661c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6629b94acceSBarry Smith 663c7afd0dbSLois Curfman McInnes Calling sequence of func: 6648d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 665c7afd0dbSLois Curfman McInnes 666313e4042SLois Curfman McInnes . f - function vector 667c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 6689b94acceSBarry Smith 6699b94acceSBarry Smith Notes: 6709b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6719b94acceSBarry Smith $ f'(x) x = -f(x), 672c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 6739b94acceSBarry Smith 67436851e7fSLois Curfman McInnes Level: beginner 67536851e7fSLois Curfman McInnes 6769b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6779b94acceSBarry Smith 678a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6799b94acceSBarry Smith @*/ 6806849ba73SBarry Smith PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 6819b94acceSBarry Smith { 6823a40ed3dSBarry Smith PetscFunctionBegin; 6834482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 6844482741eSBarry Smith PetscValidHeaderSpecific(r,VEC_COOKIE,2); 685c9780b6fSBarry Smith PetscCheckSameComm(snes,1,r,2); 686184914b5SBarry Smith 6879b94acceSBarry Smith snes->computefunction = func; 6889b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6899b94acceSBarry Smith snes->funP = ctx; 6903a40ed3dSBarry Smith PetscFunctionReturn(0); 6919b94acceSBarry Smith } 6929b94acceSBarry Smith 6933ab0aad5SBarry Smith /* --------------------------------------------------------------- */ 6943ab0aad5SBarry Smith #undef __FUNCT__ 6953ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs" 6963ab0aad5SBarry Smith /*@C 6973ab0aad5SBarry Smith SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set 6983ab0aad5SBarry Smith it assumes a zero right hand side. 6993ab0aad5SBarry Smith 7003ab0aad5SBarry Smith Collective on SNES 7013ab0aad5SBarry Smith 7023ab0aad5SBarry Smith Input Parameters: 7033ab0aad5SBarry Smith + snes - the SNES context 7043ab0aad5SBarry Smith - rhs - the right hand side vector or PETSC_NULL for a zero right hand side 7053ab0aad5SBarry Smith 7063ab0aad5SBarry Smith Level: intermediate 7073ab0aad5SBarry Smith 7083ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side 7093ab0aad5SBarry Smith 7103ab0aad5SBarry Smith .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 7113ab0aad5SBarry Smith @*/ 712dfbe8321SBarry Smith PetscErrorCode SNESSetRhs(SNES snes,Vec rhs) 7133ab0aad5SBarry Smith { 714dfbe8321SBarry Smith PetscErrorCode ierr; 7153ab0aad5SBarry Smith 7163ab0aad5SBarry Smith PetscFunctionBegin; 7173ab0aad5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7183ab0aad5SBarry Smith if (rhs) { 7193ab0aad5SBarry Smith PetscValidHeaderSpecific(rhs,VEC_COOKIE,2); 7203ab0aad5SBarry Smith PetscCheckSameComm(snes,1,rhs,2); 7213ab0aad5SBarry Smith ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr); 7223ab0aad5SBarry Smith } 7233ab0aad5SBarry Smith if (snes->afine) { 7243ab0aad5SBarry Smith ierr = VecDestroy(snes->afine);CHKERRQ(ierr); 7253ab0aad5SBarry Smith } 7263ab0aad5SBarry Smith snes->afine = rhs; 7273ab0aad5SBarry Smith PetscFunctionReturn(0); 7283ab0aad5SBarry Smith } 7293ab0aad5SBarry Smith 7304a2ae208SSatish Balay #undef __FUNCT__ 7314a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction" 7329b94acceSBarry Smith /*@ 73336851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7349b94acceSBarry Smith SNESSetFunction(). 7359b94acceSBarry Smith 736c7afd0dbSLois Curfman McInnes Collective on SNES 737c7afd0dbSLois Curfman McInnes 7389b94acceSBarry Smith Input Parameters: 739c7afd0dbSLois Curfman McInnes + snes - the SNES context 740c7afd0dbSLois Curfman McInnes - x - input vector 7419b94acceSBarry Smith 7429b94acceSBarry Smith Output Parameter: 7433638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7449b94acceSBarry Smith 7451bffabb2SLois Curfman McInnes Notes: 74636851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 74736851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 74836851e7fSLois Curfman McInnes themselves. 74936851e7fSLois Curfman McInnes 75036851e7fSLois Curfman McInnes Level: developer 75136851e7fSLois Curfman McInnes 7529b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 7539b94acceSBarry Smith 754a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 7559b94acceSBarry Smith @*/ 756dfbe8321SBarry Smith PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y) 7579b94acceSBarry Smith { 758dfbe8321SBarry Smith PetscErrorCode ierr; 7599b94acceSBarry Smith 7603a40ed3dSBarry Smith PetscFunctionBegin; 7614482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7624482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 7634482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,3); 764c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 765c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,3); 766184914b5SBarry Smith 767d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 768d64ed03dSBarry Smith PetscStackPush("SNES user function"); 7699b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 770d64ed03dSBarry Smith PetscStackPop; 7713ab0aad5SBarry Smith if (snes->afine) { 7723ab0aad5SBarry Smith PetscScalar mone = -1.0; 7733ab0aad5SBarry Smith ierr = VecAXPY(&mone,snes->afine,y);CHKERRQ(ierr); 7743ab0aad5SBarry Smith } 775ae3c334cSLois Curfman McInnes snes->nfuncs++; 776d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 7773a40ed3dSBarry Smith PetscFunctionReturn(0); 7789b94acceSBarry Smith } 7799b94acceSBarry Smith 7804a2ae208SSatish Balay #undef __FUNCT__ 7814a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian" 78262fef451SLois Curfman McInnes /*@ 78362fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 78462fef451SLois Curfman McInnes set with SNESSetJacobian(). 78562fef451SLois Curfman McInnes 786c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 787c7afd0dbSLois Curfman McInnes 78862fef451SLois Curfman McInnes Input Parameters: 789c7afd0dbSLois Curfman McInnes + snes - the SNES context 790c7afd0dbSLois Curfman McInnes - x - input vector 79162fef451SLois Curfman McInnes 79262fef451SLois Curfman McInnes Output Parameters: 793c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 79462fef451SLois Curfman McInnes . B - optional preconditioning matrix 795c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 796fee21e36SBarry Smith 79762fef451SLois Curfman McInnes Notes: 79862fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 79962fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 80062fef451SLois Curfman McInnes 80194b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 802dc5a77f8SLois Curfman McInnes flag parameter. 80362fef451SLois Curfman McInnes 80436851e7fSLois Curfman McInnes Level: developer 80536851e7fSLois Curfman McInnes 80662fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 80762fef451SLois Curfman McInnes 80894b7f48cSBarry Smith .seealso: SNESSetJacobian(), KSPSetOperators() 80962fef451SLois Curfman McInnes @*/ 810dfbe8321SBarry Smith PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8119b94acceSBarry Smith { 812dfbe8321SBarry Smith PetscErrorCode ierr; 8133a40ed3dSBarry Smith 8143a40ed3dSBarry Smith PetscFunctionBegin; 8154482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8164482741eSBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE,2); 8174482741eSBarry Smith PetscValidPointer(flg,5); 818c9780b6fSBarry Smith PetscCheckSameComm(snes,1,X,2); 8193a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 820d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 821c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 822d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 8239b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 824d64ed03dSBarry Smith PetscStackPop; 825d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 8266d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 8274482741eSBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 8284482741eSBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 8293a40ed3dSBarry Smith PetscFunctionReturn(0); 8309b94acceSBarry Smith } 8319b94acceSBarry Smith 8324a2ae208SSatish Balay #undef __FUNCT__ 8334a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian" 8349b94acceSBarry Smith /*@C 8359b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 836044dda88SLois Curfman McInnes location to store the matrix. 8379b94acceSBarry Smith 838c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 839c7afd0dbSLois Curfman McInnes 8409b94acceSBarry Smith Input Parameters: 841c7afd0dbSLois Curfman McInnes + snes - the SNES context 8429b94acceSBarry Smith . A - Jacobian matrix 8439b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 8449b94acceSBarry Smith . func - Jacobian evaluation routine 845c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 8462cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 8479b94acceSBarry Smith 8489b94acceSBarry Smith Calling sequence of func: 8498d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 8509b94acceSBarry Smith 851c7afd0dbSLois Curfman McInnes + x - input vector 8529b94acceSBarry Smith . A - Jacobian matrix 8539b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 854ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 85594b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 856c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 8579b94acceSBarry Smith 8589b94acceSBarry Smith Notes: 85994b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 8602cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 861ac21db08SLois Curfman McInnes 862ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 8639b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 8649b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 8659b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 8669b94acceSBarry Smith throughout the global iterations. 8679b94acceSBarry Smith 86836851e7fSLois Curfman McInnes Level: beginner 86936851e7fSLois Curfman McInnes 8709b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 8719b94acceSBarry Smith 87294b7f48cSBarry Smith .seealso: KSPSetOperators(), SNESSetFunction() 8739b94acceSBarry Smith @*/ 8746849ba73SBarry Smith PetscErrorCode SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 8759b94acceSBarry Smith { 876dfbe8321SBarry Smith PetscErrorCode ierr; 8773a7fca6bSBarry Smith 8783a40ed3dSBarry Smith PetscFunctionBegin; 8794482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8804482741eSBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 8814482741eSBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 882c9780b6fSBarry Smith if (A) PetscCheckSameComm(snes,1,A,2); 883c9780b6fSBarry Smith if (B) PetscCheckSameComm(snes,1,B,2); 8843a7fca6bSBarry Smith if (func) snes->computejacobian = func; 8853a7fca6bSBarry Smith if (ctx) snes->jacP = ctx; 8863a7fca6bSBarry Smith if (A) { 8873a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 8889b94acceSBarry Smith snes->jacobian = A; 8893a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 8903a7fca6bSBarry Smith } 8913a7fca6bSBarry Smith if (B) { 8923a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 8939b94acceSBarry Smith snes->jacobian_pre = B; 8943a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 8953a7fca6bSBarry Smith } 8963a40ed3dSBarry Smith PetscFunctionReturn(0); 8979b94acceSBarry Smith } 89862fef451SLois Curfman McInnes 8994a2ae208SSatish Balay #undef __FUNCT__ 9004a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian" 901c2aafc4cSSatish Balay /*@C 902b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 903b4fd4287SBarry Smith provided context for evaluating the Jacobian. 904b4fd4287SBarry Smith 905c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 906c7afd0dbSLois Curfman McInnes 907b4fd4287SBarry Smith Input Parameter: 908b4fd4287SBarry Smith . snes - the nonlinear solver context 909b4fd4287SBarry Smith 910b4fd4287SBarry Smith Output Parameters: 911c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 912b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 91300036973SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 91400036973SBarry Smith - func - location to put Jacobian function (or PETSC_NULL) 915fee21e36SBarry Smith 91636851e7fSLois Curfman McInnes Level: advanced 91736851e7fSLois Curfman McInnes 918b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 919b4fd4287SBarry Smith @*/ 9206849ba73SBarry Smith PetscErrorCode SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 921b4fd4287SBarry Smith { 9223a40ed3dSBarry Smith PetscFunctionBegin; 9234482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 924b4fd4287SBarry Smith if (A) *A = snes->jacobian; 925b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 926b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 92700036973SBarry Smith if (func) *func = snes->computejacobian; 9283a40ed3dSBarry Smith PetscFunctionReturn(0); 929b4fd4287SBarry Smith } 930b4fd4287SBarry Smith 9319b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 932dfbe8321SBarry Smith EXTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 9339b94acceSBarry Smith 9344a2ae208SSatish Balay #undef __FUNCT__ 9354a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp" 9369b94acceSBarry Smith /*@ 9379b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 938272ac6f2SLois Curfman McInnes of a nonlinear solver. 9399b94acceSBarry Smith 940fee21e36SBarry Smith Collective on SNES 941fee21e36SBarry Smith 942c7afd0dbSLois Curfman McInnes Input Parameters: 943c7afd0dbSLois Curfman McInnes + snes - the SNES context 944c7afd0dbSLois Curfman McInnes - x - the solution vector 945c7afd0dbSLois Curfman McInnes 946272ac6f2SLois Curfman McInnes Notes: 947272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 948272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 949272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 950272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 951272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 952272ac6f2SLois Curfman McInnes 95336851e7fSLois Curfman McInnes Level: advanced 95436851e7fSLois Curfman McInnes 9559b94acceSBarry Smith .keywords: SNES, nonlinear, setup 9569b94acceSBarry Smith 9579b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 9589b94acceSBarry Smith @*/ 959dfbe8321SBarry Smith PetscErrorCode SNESSetUp(SNES snes,Vec x) 9609b94acceSBarry Smith { 961dfbe8321SBarry Smith PetscErrorCode ierr; 9624b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 9633a40ed3dSBarry Smith 9643a40ed3dSBarry Smith PetscFunctionBegin; 9654482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 9664482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 967c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 9688ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 9699b94acceSBarry Smith 970b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 971c1f60f51SBarry Smith /* 972c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 973dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 974c1f60f51SBarry Smith */ 975c1f60f51SBarry Smith if (flg) { 976c1f60f51SBarry Smith Mat J; 9775a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 9785a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 9793a7fca6bSBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n"); 9803a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 9813a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 982c1f60f51SBarry Smith } 98345fc7adcSBarry Smith 984b4014bb2SMatthew Knepley #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) 98545fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 98645fc7adcSBarry Smith if (flg) { 98745fc7adcSBarry Smith Mat J; 98845fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 98945fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 99045fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 99145fc7adcSBarry Smith } 99232a4b47aSMatthew Knepley #endif 99345fc7adcSBarry Smith 994b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 995c1f60f51SBarry Smith /* 996dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 997c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 998c1f60f51SBarry Smith */ 99931615d04SLois Curfman McInnes if (flg) { 1000272ac6f2SLois Curfman McInnes Mat J; 1001b5d62d44SBarry Smith KSP ksp; 100294b7f48cSBarry Smith PC pc; 1003f3ef73ceSBarry Smith 10045a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10053a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 100693b98531SBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"); 10073a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 10083a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 10093a7fca6bSBarry Smith 1010f3ef73ceSBarry Smith /* force no preconditioner */ 101194b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1012b5d62d44SBarry Smith ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1013a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1014a9815358SBarry Smith if (!flg) { 1015f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1016272ac6f2SLois Curfman McInnes } 1017a9815358SBarry Smith } 1018f3ef73ceSBarry Smith 101929bbc08cSBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 102029bbc08cSBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 102129bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1022a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 102329bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1024a8c6a408SBarry Smith } 1025a703fe33SLois Curfman McInnes 1026a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 10274b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 10286831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 102994b7f48cSBarry Smith KSP ksp; 103094b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1031186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1032a703fe33SLois Curfman McInnes } 10334b27c08aSLois Curfman McInnes 1034a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 103582bf6240SBarry Smith snes->setupcalled = 1; 10363a40ed3dSBarry Smith PetscFunctionReturn(0); 10379b94acceSBarry Smith } 10389b94acceSBarry Smith 10394a2ae208SSatish Balay #undef __FUNCT__ 10404a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 10419b94acceSBarry Smith /*@C 10429b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 10439b94acceSBarry Smith with SNESCreate(). 10449b94acceSBarry Smith 1045c7afd0dbSLois Curfman McInnes Collective on SNES 1046c7afd0dbSLois Curfman McInnes 10479b94acceSBarry Smith Input Parameter: 10489b94acceSBarry Smith . snes - the SNES context 10499b94acceSBarry Smith 105036851e7fSLois Curfman McInnes Level: beginner 105136851e7fSLois Curfman McInnes 10529b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 10539b94acceSBarry Smith 105463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 10559b94acceSBarry Smith @*/ 1056dfbe8321SBarry Smith PetscErrorCode SNESDestroy(SNES snes) 10579b94acceSBarry Smith { 10586849ba73SBarry Smith int i; 10596849ba73SBarry Smith PetscErrorCode ierr; 10603a40ed3dSBarry Smith 10613a40ed3dSBarry Smith PetscFunctionBegin; 10624482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 10633a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1064d4bb536fSBarry Smith 1065be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 10660f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1067be0abb6dSBarry Smith 1068e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1069606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 10703a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 10713a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 10723ab0aad5SBarry Smith if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 107394b7f48cSBarry Smith ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1074522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1075b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1076b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1077b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1078b8a78c4aSBarry Smith } 1079b8a78c4aSBarry Smith } 1080b0a32e0cSBarry Smith PetscLogObjectDestroy((PetscObject)snes); 10810452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 10823a40ed3dSBarry Smith PetscFunctionReturn(0); 10839b94acceSBarry Smith } 10849b94acceSBarry Smith 10859b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 10869b94acceSBarry Smith 10874a2ae208SSatish Balay #undef __FUNCT__ 10884a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 10899b94acceSBarry Smith /*@ 1090d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 10919b94acceSBarry Smith 1092c7afd0dbSLois Curfman McInnes Collective on SNES 1093c7afd0dbSLois Curfman McInnes 10949b94acceSBarry Smith Input Parameters: 1095c7afd0dbSLois Curfman McInnes + snes - the SNES context 109633174efeSLois Curfman McInnes . atol - absolute convergence tolerance 109733174efeSLois Curfman McInnes . rtol - relative convergence tolerance 109833174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 109933174efeSLois Curfman McInnes of the change in the solution between steps 110033174efeSLois Curfman McInnes . maxit - maximum number of iterations 1101c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1102fee21e36SBarry Smith 110333174efeSLois Curfman McInnes Options Database Keys: 1104c7afd0dbSLois Curfman McInnes + -snes_atol <atol> - Sets atol 1105c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1106c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1107c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1108c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 11099b94acceSBarry Smith 1110d7a720efSLois Curfman McInnes Notes: 11119b94acceSBarry Smith The default maximum number of iterations is 50. 11129b94acceSBarry Smith The default maximum number of function evaluations is 1000. 11139b94acceSBarry Smith 111436851e7fSLois Curfman McInnes Level: intermediate 111536851e7fSLois Curfman McInnes 111633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 11179b94acceSBarry Smith 1118d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 11199b94acceSBarry Smith @*/ 1120dfbe8321SBarry Smith PetscErrorCode SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf) 11219b94acceSBarry Smith { 11223a40ed3dSBarry Smith PetscFunctionBegin; 11234482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1124d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1125d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1126d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1127d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1128d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11293a40ed3dSBarry Smith PetscFunctionReturn(0); 11309b94acceSBarry Smith } 11319b94acceSBarry Smith 11324a2ae208SSatish Balay #undef __FUNCT__ 11334a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11349b94acceSBarry Smith /*@ 113533174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 113633174efeSLois Curfman McInnes 1137c7afd0dbSLois Curfman McInnes Not Collective 1138c7afd0dbSLois Curfman McInnes 113933174efeSLois Curfman McInnes Input Parameters: 1140c7afd0dbSLois Curfman McInnes + snes - the SNES context 114133174efeSLois Curfman McInnes . atol - absolute convergence tolerance 114233174efeSLois Curfman McInnes . rtol - relative convergence tolerance 114333174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 114433174efeSLois Curfman McInnes of the change in the solution between steps 114533174efeSLois Curfman McInnes . maxit - maximum number of iterations 1146c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1147fee21e36SBarry Smith 114833174efeSLois Curfman McInnes Notes: 114933174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 115033174efeSLois Curfman McInnes 115136851e7fSLois Curfman McInnes Level: intermediate 115236851e7fSLois Curfman McInnes 115333174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 115433174efeSLois Curfman McInnes 115533174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 115633174efeSLois Curfman McInnes @*/ 1157dfbe8321SBarry Smith PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf) 115833174efeSLois Curfman McInnes { 11593a40ed3dSBarry Smith PetscFunctionBegin; 11604482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 116133174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 116233174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 116333174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 116433174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 116533174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 11663a40ed3dSBarry Smith PetscFunctionReturn(0); 116733174efeSLois Curfman McInnes } 116833174efeSLois Curfman McInnes 11694a2ae208SSatish Balay #undef __FUNCT__ 11704a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 117133174efeSLois Curfman McInnes /*@ 11729b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 11739b94acceSBarry Smith 1174fee21e36SBarry Smith Collective on SNES 1175fee21e36SBarry Smith 1176c7afd0dbSLois Curfman McInnes Input Parameters: 1177c7afd0dbSLois Curfman McInnes + snes - the SNES context 1178c7afd0dbSLois Curfman McInnes - tol - tolerance 1179c7afd0dbSLois Curfman McInnes 11809b94acceSBarry Smith Options Database Key: 1181c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 11829b94acceSBarry Smith 118336851e7fSLois Curfman McInnes Level: intermediate 118436851e7fSLois Curfman McInnes 11859b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 11869b94acceSBarry Smith 1187d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 11889b94acceSBarry Smith @*/ 1189dfbe8321SBarry Smith PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 11909b94acceSBarry Smith { 11913a40ed3dSBarry Smith PetscFunctionBegin; 11924482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 11939b94acceSBarry Smith snes->deltatol = tol; 11943a40ed3dSBarry Smith PetscFunctionReturn(0); 11959b94acceSBarry Smith } 11969b94acceSBarry Smith 1197df9fa365SBarry Smith /* 1198df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1199df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1200df9fa365SBarry Smith macros instead of functions 1201df9fa365SBarry Smith */ 12024a2ae208SSatish Balay #undef __FUNCT__ 12034a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 1204dfbe8321SBarry Smith PetscErrorCode SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx) 1205ce1608b8SBarry Smith { 1206dfbe8321SBarry Smith PetscErrorCode ierr; 1207ce1608b8SBarry Smith 1208ce1608b8SBarry Smith PetscFunctionBegin; 12094482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1210ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1211ce1608b8SBarry Smith PetscFunctionReturn(0); 1212ce1608b8SBarry Smith } 1213ce1608b8SBarry Smith 12144a2ae208SSatish Balay #undef __FUNCT__ 12154a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 1216dfbe8321SBarry Smith PetscErrorCode SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1217df9fa365SBarry Smith { 1218dfbe8321SBarry Smith PetscErrorCode ierr; 1219df9fa365SBarry Smith 1220df9fa365SBarry Smith PetscFunctionBegin; 1221df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1222df9fa365SBarry Smith PetscFunctionReturn(0); 1223df9fa365SBarry Smith } 1224df9fa365SBarry Smith 12254a2ae208SSatish Balay #undef __FUNCT__ 12264a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 1227dfbe8321SBarry Smith PetscErrorCode SNESLGMonitorDestroy(PetscDrawLG draw) 1228df9fa365SBarry Smith { 1229dfbe8321SBarry Smith PetscErrorCode ierr; 1230df9fa365SBarry Smith 1231df9fa365SBarry Smith PetscFunctionBegin; 1232df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1233df9fa365SBarry Smith PetscFunctionReturn(0); 1234df9fa365SBarry Smith } 1235df9fa365SBarry Smith 12369b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12379b94acceSBarry Smith 12384a2ae208SSatish Balay #undef __FUNCT__ 12394a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12409b94acceSBarry Smith /*@C 12415cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12429b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12439b94acceSBarry Smith progress. 12449b94acceSBarry Smith 1245fee21e36SBarry Smith Collective on SNES 1246fee21e36SBarry Smith 1247c7afd0dbSLois Curfman McInnes Input Parameters: 1248c7afd0dbSLois Curfman McInnes + snes - the SNES context 1249c7afd0dbSLois Curfman McInnes . func - monitoring routine 1250b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1251b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desitre) 1252b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1253b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 12549b94acceSBarry Smith 1255c7afd0dbSLois Curfman McInnes Calling sequence of func: 1256329f5518SBarry Smith $ int func(SNES snes,int its, PetscReal norm,void *mctx) 1257c7afd0dbSLois Curfman McInnes 1258c7afd0dbSLois Curfman McInnes + snes - the SNES context 1259c7afd0dbSLois Curfman McInnes . its - iteration number 1260c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 126140a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 12629b94acceSBarry Smith 12639665c990SLois Curfman McInnes Options Database Keys: 1264c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1265c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1266c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1267c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1268c7afd0dbSLois Curfman McInnes been hardwired into a code by 1269c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1270c7afd0dbSLois Curfman McInnes does not cancel those set via 1271c7afd0dbSLois Curfman McInnes the options database. 12729665c990SLois Curfman McInnes 1273639f9d9dSBarry Smith Notes: 12746bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 12756bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 12766bc08f3fSLois Curfman McInnes order in which they were set. 1277639f9d9dSBarry Smith 127836851e7fSLois Curfman McInnes Level: intermediate 127936851e7fSLois Curfman McInnes 12809b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 12819b94acceSBarry Smith 12825cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 12839b94acceSBarry Smith @*/ 12846849ba73SBarry Smith PetscErrorCode SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,int,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 12859b94acceSBarry Smith { 12863a40ed3dSBarry Smith PetscFunctionBegin; 12874482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1288639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 128929bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1290639f9d9dSBarry Smith } 1291639f9d9dSBarry Smith 1292639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1293b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1294639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 12953a40ed3dSBarry Smith PetscFunctionReturn(0); 12969b94acceSBarry Smith } 12979b94acceSBarry Smith 12984a2ae208SSatish Balay #undef __FUNCT__ 12994a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 13005cd90555SBarry Smith /*@C 13015cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 13025cd90555SBarry Smith 1303c7afd0dbSLois Curfman McInnes Collective on SNES 1304c7afd0dbSLois Curfman McInnes 13055cd90555SBarry Smith Input Parameters: 13065cd90555SBarry Smith . snes - the SNES context 13075cd90555SBarry Smith 13081a480d89SAdministrator Options Database Key: 1309c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1310c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1311c7afd0dbSLois Curfman McInnes set via the options database 13125cd90555SBarry Smith 13135cd90555SBarry Smith Notes: 13145cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 13155cd90555SBarry Smith 131636851e7fSLois Curfman McInnes Level: intermediate 131736851e7fSLois Curfman McInnes 13185cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 13195cd90555SBarry Smith 13205cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 13215cd90555SBarry Smith @*/ 1322dfbe8321SBarry Smith PetscErrorCode SNESClearMonitor(SNES snes) 13235cd90555SBarry Smith { 13245cd90555SBarry Smith PetscFunctionBegin; 13254482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13265cd90555SBarry Smith snes->numbermonitors = 0; 13275cd90555SBarry Smith PetscFunctionReturn(0); 13285cd90555SBarry Smith } 13295cd90555SBarry Smith 13304a2ae208SSatish Balay #undef __FUNCT__ 13314a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13329b94acceSBarry Smith /*@C 13339b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13349b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13359b94acceSBarry Smith 1336fee21e36SBarry Smith Collective on SNES 1337fee21e36SBarry Smith 1338c7afd0dbSLois Curfman McInnes Input Parameters: 1339c7afd0dbSLois Curfman McInnes + snes - the SNES context 1340c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1341c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1342c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 13439b94acceSBarry Smith 1344c7afd0dbSLois Curfman McInnes Calling sequence of func: 13456849ba73SBarry Smith $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1346c7afd0dbSLois Curfman McInnes 1347c7afd0dbSLois Curfman McInnes + snes - the SNES context 1348c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1349184914b5SBarry Smith . reason - reason for convergence/divergence 1350c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 13514b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 13524b27c08aSLois Curfman McInnes - f - 2-norm of function 13539b94acceSBarry Smith 135436851e7fSLois Curfman McInnes Level: advanced 135536851e7fSLois Curfman McInnes 13569b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13579b94acceSBarry Smith 13584b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 13599b94acceSBarry Smith @*/ 13606849ba73SBarry Smith PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 13619b94acceSBarry Smith { 13623a40ed3dSBarry Smith PetscFunctionBegin; 13634482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13649b94acceSBarry Smith (snes)->converged = func; 13659b94acceSBarry Smith (snes)->cnvP = cctx; 13663a40ed3dSBarry Smith PetscFunctionReturn(0); 13679b94acceSBarry Smith } 13689b94acceSBarry Smith 13694a2ae208SSatish Balay #undef __FUNCT__ 13704a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 1371184914b5SBarry Smith /*@C 1372184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1373184914b5SBarry Smith 1374184914b5SBarry Smith Not Collective 1375184914b5SBarry Smith 1376184914b5SBarry Smith Input Parameter: 1377184914b5SBarry Smith . snes - the SNES context 1378184914b5SBarry Smith 1379184914b5SBarry Smith Output Parameter: 1380e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1381184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1382184914b5SBarry Smith 1383184914b5SBarry Smith Level: intermediate 1384184914b5SBarry Smith 1385184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1386184914b5SBarry Smith 1387184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1388184914b5SBarry Smith 13894b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1390184914b5SBarry Smith @*/ 1391dfbe8321SBarry Smith PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1392184914b5SBarry Smith { 1393184914b5SBarry Smith PetscFunctionBegin; 13944482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13954482741eSBarry Smith PetscValidPointer(reason,2); 1396184914b5SBarry Smith *reason = snes->reason; 1397184914b5SBarry Smith PetscFunctionReturn(0); 1398184914b5SBarry Smith } 1399184914b5SBarry Smith 14004a2ae208SSatish Balay #undef __FUNCT__ 14014a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1402c9005455SLois Curfman McInnes /*@ 1403c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1404c9005455SLois Curfman McInnes 1405fee21e36SBarry Smith Collective on SNES 1406fee21e36SBarry Smith 1407c7afd0dbSLois Curfman McInnes Input Parameters: 1408c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1409c7afd0dbSLois Curfman McInnes . a - array to hold history 1410cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1411758f92a0SBarry Smith . na - size of a and its 141264731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1413758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1414c7afd0dbSLois Curfman McInnes 1415c9005455SLois Curfman McInnes Notes: 14164b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1417c9005455SLois Curfman McInnes at each step. 1418c9005455SLois Curfman McInnes 1419c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1420c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1421c9005455SLois Curfman McInnes during the section of code that is being timed. 1422c9005455SLois Curfman McInnes 142336851e7fSLois Curfman McInnes Level: intermediate 142436851e7fSLois Curfman McInnes 1425c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1426758f92a0SBarry Smith 142708405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1428758f92a0SBarry Smith 1429c9005455SLois Curfman McInnes @*/ 1430dfbe8321SBarry Smith PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],int *its,int na,PetscTruth reset) 1431c9005455SLois Curfman McInnes { 14323a40ed3dSBarry Smith PetscFunctionBegin; 14334482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14344482741eSBarry Smith if (na) PetscValidScalarPointer(a,2); 1435c9005455SLois Curfman McInnes snes->conv_hist = a; 1436758f92a0SBarry Smith snes->conv_hist_its = its; 1437758f92a0SBarry Smith snes->conv_hist_max = na; 1438758f92a0SBarry Smith snes->conv_hist_reset = reset; 1439758f92a0SBarry Smith PetscFunctionReturn(0); 1440758f92a0SBarry Smith } 1441758f92a0SBarry Smith 14424a2ae208SSatish Balay #undef __FUNCT__ 14434a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 14440c4c9dddSBarry Smith /*@C 1445758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1446758f92a0SBarry Smith 1447758f92a0SBarry Smith Collective on SNES 1448758f92a0SBarry Smith 1449758f92a0SBarry Smith Input Parameter: 1450758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1451758f92a0SBarry Smith 1452758f92a0SBarry Smith Output Parameters: 1453758f92a0SBarry Smith . a - array to hold history 1454758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1455758f92a0SBarry Smith negative if not converged) for each solve. 1456758f92a0SBarry Smith - na - size of a and its 1457758f92a0SBarry Smith 1458758f92a0SBarry Smith Notes: 1459758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1460758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1461758f92a0SBarry Smith 1462758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1463758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1464758f92a0SBarry Smith during the section of code that is being timed. 1465758f92a0SBarry Smith 1466758f92a0SBarry Smith Level: intermediate 1467758f92a0SBarry Smith 1468758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1469758f92a0SBarry Smith 1470758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1471758f92a0SBarry Smith 1472758f92a0SBarry Smith @*/ 1473dfbe8321SBarry Smith PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],int *its[],int *na) 1474758f92a0SBarry Smith { 1475758f92a0SBarry Smith PetscFunctionBegin; 14764482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1477758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1478758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1479758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 14803a40ed3dSBarry Smith PetscFunctionReturn(0); 1481c9005455SLois Curfman McInnes } 1482c9005455SLois Curfman McInnes 1483e74ef692SMatthew Knepley #undef __FUNCT__ 1484e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC" 1485ac226902SBarry Smith /*@C 148676b2cf59SMatthew Knepley SNESSetRhsBC - Sets the function which applies boundary conditions 148776b2cf59SMatthew Knepley to the Rhs of each system. 148876b2cf59SMatthew Knepley 148976b2cf59SMatthew Knepley Collective on SNES 149076b2cf59SMatthew Knepley 149176b2cf59SMatthew Knepley Input Parameters: 149276b2cf59SMatthew Knepley . snes - The nonlinear solver context 149376b2cf59SMatthew Knepley . func - The function 149476b2cf59SMatthew Knepley 149576b2cf59SMatthew Knepley Calling sequence of func: 149676b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx); 149776b2cf59SMatthew Knepley 149876b2cf59SMatthew Knepley . rhs - The current rhs vector 149976b2cf59SMatthew Knepley . ctx - The user-context 150076b2cf59SMatthew Knepley 150176b2cf59SMatthew Knepley Level: intermediate 150276b2cf59SMatthew Knepley 150376b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 150476b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate() 150576b2cf59SMatthew Knepley @*/ 15066849ba73SBarry Smith PetscErrorCode SNESSetRhsBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 150776b2cf59SMatthew Knepley { 150876b2cf59SMatthew Knepley PetscFunctionBegin; 15094482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 151076b2cf59SMatthew Knepley snes->applyrhsbc = func; 151176b2cf59SMatthew Knepley PetscFunctionReturn(0); 151276b2cf59SMatthew Knepley } 151376b2cf59SMatthew Knepley 1514e74ef692SMatthew Knepley #undef __FUNCT__ 1515e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC" 151676b2cf59SMatthew Knepley /*@ 151776b2cf59SMatthew Knepley SNESDefaultRhsBC - The default boundary condition function which does nothing. 151876b2cf59SMatthew Knepley 151976b2cf59SMatthew Knepley Not collective 152076b2cf59SMatthew Knepley 152176b2cf59SMatthew Knepley Input Parameters: 152276b2cf59SMatthew Knepley . snes - The nonlinear solver context 152376b2cf59SMatthew Knepley . rhs - The Rhs 152476b2cf59SMatthew Knepley . ctx - The user-context 152576b2cf59SMatthew Knepley 152676b2cf59SMatthew Knepley Level: beginner 152776b2cf59SMatthew Knepley 152876b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 152976b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate() 153076b2cf59SMatthew Knepley @*/ 1531dfbe8321SBarry Smith PetscErrorCode SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx) 153276b2cf59SMatthew Knepley { 153376b2cf59SMatthew Knepley PetscFunctionBegin; 153476b2cf59SMatthew Knepley PetscFunctionReturn(0); 153576b2cf59SMatthew Knepley } 153676b2cf59SMatthew Knepley 1537e74ef692SMatthew Knepley #undef __FUNCT__ 1538e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC" 1539ac226902SBarry Smith /*@C 154076b2cf59SMatthew Knepley SNESSetSolutionBC - Sets the function which applies boundary conditions 154176b2cf59SMatthew Knepley to the solution of each system. 154276b2cf59SMatthew Knepley 154376b2cf59SMatthew Knepley Collective on SNES 154476b2cf59SMatthew Knepley 154576b2cf59SMatthew Knepley Input Parameters: 154676b2cf59SMatthew Knepley . snes - The nonlinear solver context 154776b2cf59SMatthew Knepley . func - The function 154876b2cf59SMatthew Knepley 154976b2cf59SMatthew Knepley Calling sequence of func: 15509d06fb6bSMatthew Knepley . func (SNES snes, Vec rsol, void *ctx); 155176b2cf59SMatthew Knepley 155276b2cf59SMatthew Knepley . sol - The current solution vector 155376b2cf59SMatthew Knepley . ctx - The user-context 155476b2cf59SMatthew Knepley 155576b2cf59SMatthew Knepley Level: intermediate 155676b2cf59SMatthew Knepley 155776b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 155876b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate() 155976b2cf59SMatthew Knepley @*/ 15606849ba73SBarry Smith PetscErrorCode SNESSetSolutionBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 156176b2cf59SMatthew Knepley { 156276b2cf59SMatthew Knepley PetscFunctionBegin; 15634482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 156476b2cf59SMatthew Knepley snes->applysolbc = func; 156576b2cf59SMatthew Knepley PetscFunctionReturn(0); 156676b2cf59SMatthew Knepley } 156776b2cf59SMatthew Knepley 1568e74ef692SMatthew Knepley #undef __FUNCT__ 1569e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC" 157076b2cf59SMatthew Knepley /*@ 157176b2cf59SMatthew Knepley SNESDefaultSolutionBC - The default boundary condition function which does nothing. 157276b2cf59SMatthew Knepley 157376b2cf59SMatthew Knepley Not collective 157476b2cf59SMatthew Knepley 157576b2cf59SMatthew Knepley Input Parameters: 157676b2cf59SMatthew Knepley . snes - The nonlinear solver context 157776b2cf59SMatthew Knepley . sol - The solution 157876b2cf59SMatthew Knepley . ctx - The user-context 157976b2cf59SMatthew Knepley 158076b2cf59SMatthew Knepley Level: beginner 158176b2cf59SMatthew Knepley 158276b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 158376b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate() 158476b2cf59SMatthew Knepley @*/ 1585dfbe8321SBarry Smith PetscErrorCode SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx) 158676b2cf59SMatthew Knepley { 158776b2cf59SMatthew Knepley PetscFunctionBegin; 158876b2cf59SMatthew Knepley PetscFunctionReturn(0); 158976b2cf59SMatthew Knepley } 159076b2cf59SMatthew Knepley 1591e74ef692SMatthew Knepley #undef __FUNCT__ 1592e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 1593ac226902SBarry Smith /*@C 159476b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 159576b2cf59SMatthew Knepley at the beginning of every step of the iteration. 159676b2cf59SMatthew Knepley 159776b2cf59SMatthew Knepley Collective on SNES 159876b2cf59SMatthew Knepley 159976b2cf59SMatthew Knepley Input Parameters: 160076b2cf59SMatthew Knepley . snes - The nonlinear solver context 160176b2cf59SMatthew Knepley . func - The function 160276b2cf59SMatthew Knepley 160376b2cf59SMatthew Knepley Calling sequence of func: 160476b2cf59SMatthew Knepley . func (TS ts, int step); 160576b2cf59SMatthew Knepley 160676b2cf59SMatthew Knepley . step - The current step of the iteration 160776b2cf59SMatthew Knepley 160876b2cf59SMatthew Knepley Level: intermediate 160976b2cf59SMatthew Knepley 161076b2cf59SMatthew Knepley .keywords: SNES, update 161176b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 161276b2cf59SMatthew Knepley @*/ 16136849ba73SBarry Smith PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, int)) 161476b2cf59SMatthew Knepley { 161576b2cf59SMatthew Knepley PetscFunctionBegin; 16164482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 161776b2cf59SMatthew Knepley snes->update = func; 161876b2cf59SMatthew Knepley PetscFunctionReturn(0); 161976b2cf59SMatthew Knepley } 162076b2cf59SMatthew Knepley 1621e74ef692SMatthew Knepley #undef __FUNCT__ 1622e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 162376b2cf59SMatthew Knepley /*@ 162476b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 162576b2cf59SMatthew Knepley 162676b2cf59SMatthew Knepley Not collective 162776b2cf59SMatthew Knepley 162876b2cf59SMatthew Knepley Input Parameters: 162976b2cf59SMatthew Knepley . snes - The nonlinear solver context 163076b2cf59SMatthew Knepley . step - The current step of the iteration 163176b2cf59SMatthew Knepley 1632205452f4SMatthew Knepley Level: intermediate 1633205452f4SMatthew Knepley 163476b2cf59SMatthew Knepley .keywords: SNES, update 163576b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 163676b2cf59SMatthew Knepley @*/ 1637dfbe8321SBarry Smith PetscErrorCode SNESDefaultUpdate(SNES snes, int step) 163876b2cf59SMatthew Knepley { 163976b2cf59SMatthew Knepley PetscFunctionBegin; 164076b2cf59SMatthew Knepley PetscFunctionReturn(0); 164176b2cf59SMatthew Knepley } 164276b2cf59SMatthew Knepley 16434a2ae208SSatish Balay #undef __FUNCT__ 16444a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 16459b94acceSBarry Smith /* 16469b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 16479b94acceSBarry Smith positive parameter delta. 16489b94acceSBarry Smith 16499b94acceSBarry Smith Input Parameters: 1650c7afd0dbSLois Curfman McInnes + snes - the SNES context 16519b94acceSBarry Smith . y - approximate solution of linear system 16529b94acceSBarry Smith . fnorm - 2-norm of current function 1653c7afd0dbSLois Curfman McInnes - delta - trust region size 16549b94acceSBarry Smith 16559b94acceSBarry Smith Output Parameters: 1656c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 16579b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 16589b94acceSBarry Smith region, and exceeds zero otherwise. 1659c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 16609b94acceSBarry Smith 16619b94acceSBarry Smith Note: 16624b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 16639b94acceSBarry Smith is set to be the maximum allowable step size. 16649b94acceSBarry Smith 16659b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16669b94acceSBarry Smith */ 1667dfbe8321SBarry Smith PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 16689b94acceSBarry Smith { 1669064f8208SBarry Smith PetscReal nrm; 1670ea709b57SSatish Balay PetscScalar cnorm; 1671dfbe8321SBarry Smith PetscErrorCode ierr; 16723a40ed3dSBarry Smith 16733a40ed3dSBarry Smith PetscFunctionBegin; 16744482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16754482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1676c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,2); 1677184914b5SBarry Smith 1678064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1679064f8208SBarry Smith if (nrm > *delta) { 1680064f8208SBarry Smith nrm = *delta/nrm; 1681064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1682064f8208SBarry Smith cnorm = nrm; 1683329f5518SBarry Smith ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 16849b94acceSBarry Smith *ynorm = *delta; 16859b94acceSBarry Smith } else { 16869b94acceSBarry Smith *gpnorm = 0.0; 1687064f8208SBarry Smith *ynorm = nrm; 16889b94acceSBarry Smith } 16893a40ed3dSBarry Smith PetscFunctionReturn(0); 16909b94acceSBarry Smith } 16919b94acceSBarry Smith 16925968eb51SBarry Smith static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero", 16935968eb51SBarry Smith "not currently used", 16945968eb51SBarry Smith "line search failed", 16955968eb51SBarry Smith "reach maximum number of iterations", 16965968eb51SBarry Smith "function norm became NaN (not a number)", 16975968eb51SBarry Smith "not currently used", 16985968eb51SBarry Smith "number of function computations exceeded", 16995968eb51SBarry Smith "not currently used", 17005968eb51SBarry Smith "still iterating", 17015968eb51SBarry Smith "not currently used", 17025968eb51SBarry Smith "absolute size of function norm", 17035968eb51SBarry Smith "relative decrease in function norm", 17045968eb51SBarry Smith "step size is small", 17055968eb51SBarry Smith "not currently used", 17065968eb51SBarry Smith "not currently used", 17075968eb51SBarry Smith "small trust region"}; 17085968eb51SBarry Smith 17094a2ae208SSatish Balay #undef __FUNCT__ 17104a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 17119b94acceSBarry Smith /*@ 17129b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 171363a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 17149b94acceSBarry Smith 1715c7afd0dbSLois Curfman McInnes Collective on SNES 1716c7afd0dbSLois Curfman McInnes 1717b2002411SLois Curfman McInnes Input Parameters: 1718c7afd0dbSLois Curfman McInnes + snes - the SNES context 1719c7afd0dbSLois Curfman McInnes - x - the solution vector 17209b94acceSBarry Smith 1721b2002411SLois Curfman McInnes Notes: 17228ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 17238ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 17248ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 17258ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 17268ddd3da0SLois Curfman McInnes 172736851e7fSLois Curfman McInnes Level: beginner 172836851e7fSLois Curfman McInnes 17299b94acceSBarry Smith .keywords: SNES, nonlinear, solve 17309b94acceSBarry Smith 17313ab0aad5SBarry Smith .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs() 17329b94acceSBarry Smith @*/ 1733dfbe8321SBarry Smith PetscErrorCode SNESSolve(SNES snes,Vec x) 17349b94acceSBarry Smith { 1735dfbe8321SBarry Smith PetscErrorCode ierr; 1736f1af5d2fSBarry Smith PetscTruth flg; 1737052efed2SBarry Smith 17383a40ed3dSBarry Smith PetscFunctionBegin; 17394482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17404482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1741c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 1742*1302d50aSBarry Smith if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 174374637425SBarry Smith 174482bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1745c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 1746758f92a0SBarry Smith if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 1747d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 174850ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1749c9780b6fSBarry Smith ierr = (*(snes)->solve)(snes);CHKERRQ(ierr); 1750d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1751b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 17528b179ff8SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1753da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1754da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 17555968eb51SBarry Smith if (snes->printreason) { 17565968eb51SBarry Smith if (snes->reason > 0) { 17575968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 17585968eb51SBarry Smith } else { 17595968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 17605968eb51SBarry Smith } 17615968eb51SBarry Smith } 17625968eb51SBarry Smith 17633a40ed3dSBarry Smith PetscFunctionReturn(0); 17649b94acceSBarry Smith } 17659b94acceSBarry Smith 17669b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 17679b94acceSBarry Smith 17684a2ae208SSatish Balay #undef __FUNCT__ 17694a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 177082bf6240SBarry Smith /*@C 17714b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 17729b94acceSBarry Smith 1773fee21e36SBarry Smith Collective on SNES 1774fee21e36SBarry Smith 1775c7afd0dbSLois Curfman McInnes Input Parameters: 1776c7afd0dbSLois Curfman McInnes + snes - the SNES context 1777454a90a3SBarry Smith - type - a known method 1778c7afd0dbSLois Curfman McInnes 1779c7afd0dbSLois Curfman McInnes Options Database Key: 1780454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1781c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1782ae12b187SLois Curfman McInnes 17839b94acceSBarry Smith Notes: 1784e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 17854b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1786c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17874b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1788c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17899b94acceSBarry Smith 1790ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1791ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1792ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1793ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1794ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1795ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1796ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1797ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1798ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1799b0a32e0cSBarry Smith appropriate method. 180036851e7fSLois Curfman McInnes 180136851e7fSLois Curfman McInnes Level: intermediate 1802a703fe33SLois Curfman McInnes 1803454a90a3SBarry Smith .keywords: SNES, set, type 1804435da068SBarry Smith 1805435da068SBarry Smith .seealso: SNESType, SNESCreate() 1806435da068SBarry Smith 18079b94acceSBarry Smith @*/ 1808dfbe8321SBarry Smith PetscErrorCode SNESSetType(SNES snes,const SNESType type) 18099b94acceSBarry Smith { 1810dfbe8321SBarry Smith PetscErrorCode ierr,(*r)(SNES); 18116831982aSBarry Smith PetscTruth match; 18123a40ed3dSBarry Smith 18133a40ed3dSBarry Smith PetscFunctionBegin; 18144482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18154482741eSBarry Smith PetscValidCharPointer(type,2); 181682bf6240SBarry Smith 18176831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 18180f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 181982bf6240SBarry Smith 182082bf6240SBarry Smith if (snes->setupcalled) { 1821e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 182282bf6240SBarry Smith snes->data = 0; 182382bf6240SBarry Smith } 18247d1a2b2bSBarry Smith 18259b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 182682bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1827b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1828958c9bccSBarry Smith if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1829606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 183082bf6240SBarry Smith snes->data = 0; 18313a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1832454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 18333a40ed3dSBarry Smith PetscFunctionReturn(0); 18349b94acceSBarry Smith } 18359b94acceSBarry Smith 1836a847f771SSatish Balay 18379b94acceSBarry Smith /* --------------------------------------------------------------------- */ 18384a2ae208SSatish Balay #undef __FUNCT__ 18394a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 18409b94acceSBarry Smith /*@C 18419b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1842f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 18439b94acceSBarry Smith 1844fee21e36SBarry Smith Not Collective 1845fee21e36SBarry Smith 184636851e7fSLois Curfman McInnes Level: advanced 184736851e7fSLois Curfman McInnes 18489b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 18499b94acceSBarry Smith 18509b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 18519b94acceSBarry Smith @*/ 1852dfbe8321SBarry Smith PetscErrorCode SNESRegisterDestroy(void) 18539b94acceSBarry Smith { 1854dfbe8321SBarry Smith PetscErrorCode ierr; 185582bf6240SBarry Smith 18563a40ed3dSBarry Smith PetscFunctionBegin; 185782bf6240SBarry Smith if (SNESList) { 1858b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 185982bf6240SBarry Smith SNESList = 0; 18609b94acceSBarry Smith } 18614c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 18623a40ed3dSBarry Smith PetscFunctionReturn(0); 18639b94acceSBarry Smith } 18649b94acceSBarry Smith 18654a2ae208SSatish Balay #undef __FUNCT__ 18664a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 18679b94acceSBarry Smith /*@C 18689a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 18699b94acceSBarry Smith 1870c7afd0dbSLois Curfman McInnes Not Collective 1871c7afd0dbSLois Curfman McInnes 18729b94acceSBarry Smith Input Parameter: 18734b0e389bSBarry Smith . snes - nonlinear solver context 18749b94acceSBarry Smith 18759b94acceSBarry Smith Output Parameter: 18763a7fca6bSBarry Smith . type - SNES method (a character string) 18779b94acceSBarry Smith 187836851e7fSLois Curfman McInnes Level: intermediate 187936851e7fSLois Curfman McInnes 1880454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 18819b94acceSBarry Smith @*/ 1882dfbe8321SBarry Smith PetscErrorCode SNESGetType(SNES snes,SNESType *type) 18839b94acceSBarry Smith { 18843a40ed3dSBarry Smith PetscFunctionBegin; 18854482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18864482741eSBarry Smith PetscValidPointer(type,2); 1887454a90a3SBarry Smith *type = snes->type_name; 18883a40ed3dSBarry Smith PetscFunctionReturn(0); 18899b94acceSBarry Smith } 18909b94acceSBarry Smith 18914a2ae208SSatish Balay #undef __FUNCT__ 18924a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 18939b94acceSBarry Smith /*@C 18949b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 18959b94acceSBarry Smith stored. 18969b94acceSBarry Smith 1897c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1898c7afd0dbSLois Curfman McInnes 18999b94acceSBarry Smith Input Parameter: 19009b94acceSBarry Smith . snes - the SNES context 19019b94acceSBarry Smith 19029b94acceSBarry Smith Output Parameter: 19039b94acceSBarry Smith . x - the solution 19049b94acceSBarry Smith 190536851e7fSLois Curfman McInnes Level: advanced 190636851e7fSLois Curfman McInnes 19079b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 19089b94acceSBarry Smith 19094b27c08aSLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetSolutionUpdate() 19109b94acceSBarry Smith @*/ 1911dfbe8321SBarry Smith PetscErrorCode SNESGetSolution(SNES snes,Vec *x) 19129b94acceSBarry Smith { 19133a40ed3dSBarry Smith PetscFunctionBegin; 19144482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19154482741eSBarry Smith PetscValidPointer(x,2); 19169b94acceSBarry Smith *x = snes->vec_sol_always; 19173a40ed3dSBarry Smith PetscFunctionReturn(0); 19189b94acceSBarry Smith } 19199b94acceSBarry Smith 19204a2ae208SSatish Balay #undef __FUNCT__ 19214a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 19229b94acceSBarry Smith /*@C 19239b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 19249b94acceSBarry Smith stored. 19259b94acceSBarry Smith 1926c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1927c7afd0dbSLois Curfman McInnes 19289b94acceSBarry Smith Input Parameter: 19299b94acceSBarry Smith . snes - the SNES context 19309b94acceSBarry Smith 19319b94acceSBarry Smith Output Parameter: 19329b94acceSBarry Smith . x - the solution update 19339b94acceSBarry Smith 193436851e7fSLois Curfman McInnes Level: advanced 193536851e7fSLois Curfman McInnes 19369b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 19379b94acceSBarry Smith 19389b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 19399b94acceSBarry Smith @*/ 1940dfbe8321SBarry Smith PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x) 19419b94acceSBarry Smith { 19423a40ed3dSBarry Smith PetscFunctionBegin; 19434482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19444482741eSBarry Smith PetscValidPointer(x,2); 19459b94acceSBarry Smith *x = snes->vec_sol_update_always; 19463a40ed3dSBarry Smith PetscFunctionReturn(0); 19479b94acceSBarry Smith } 19489b94acceSBarry Smith 19494a2ae208SSatish Balay #undef __FUNCT__ 19504a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 19519b94acceSBarry Smith /*@C 19523638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 19539b94acceSBarry Smith 1954c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1955c7afd0dbSLois Curfman McInnes 19569b94acceSBarry Smith Input Parameter: 19579b94acceSBarry Smith . snes - the SNES context 19589b94acceSBarry Smith 19599b94acceSBarry Smith Output Parameter: 19607bf4e008SBarry Smith + r - the function (or PETSC_NULL) 196100036973SBarry Smith . ctx - the function context (or PETSC_NULL) 196200036973SBarry Smith - func - the function (or PETSC_NULL) 19639b94acceSBarry Smith 196436851e7fSLois Curfman McInnes Level: advanced 196536851e7fSLois Curfman McInnes 1966a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 19679b94acceSBarry Smith 19684b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 19699b94acceSBarry Smith @*/ 19706849ba73SBarry Smith PetscErrorCode SNESGetFunction(SNES snes,Vec *r,void **ctx,PetscErrorCode (**func)(SNES,Vec,Vec,void*)) 19719b94acceSBarry Smith { 19723a40ed3dSBarry Smith PetscFunctionBegin; 19734482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19747bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 19757bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 197600036973SBarry Smith if (func) *func = snes->computefunction; 19773a40ed3dSBarry Smith PetscFunctionReturn(0); 19789b94acceSBarry Smith } 19799b94acceSBarry Smith 19804a2ae208SSatish Balay #undef __FUNCT__ 19814a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 19823c7409f5SSatish Balay /*@C 19833c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1984d850072dSLois Curfman McInnes SNES options in the database. 19853c7409f5SSatish Balay 1986fee21e36SBarry Smith Collective on SNES 1987fee21e36SBarry Smith 1988c7afd0dbSLois Curfman McInnes Input Parameter: 1989c7afd0dbSLois Curfman McInnes + snes - the SNES context 1990c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1991c7afd0dbSLois Curfman McInnes 1992d850072dSLois Curfman McInnes Notes: 1993a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1994c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1995d850072dSLois Curfman McInnes 199636851e7fSLois Curfman McInnes Level: advanced 199736851e7fSLois Curfman McInnes 19983c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1999a86d99e1SLois Curfman McInnes 2000a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 20013c7409f5SSatish Balay @*/ 2002dfbe8321SBarry Smith PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[]) 20033c7409f5SSatish Balay { 2004dfbe8321SBarry Smith PetscErrorCode ierr; 20053c7409f5SSatish Balay 20063a40ed3dSBarry Smith PetscFunctionBegin; 20074482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2008639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 200994b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20103a40ed3dSBarry Smith PetscFunctionReturn(0); 20113c7409f5SSatish Balay } 20123c7409f5SSatish Balay 20134a2ae208SSatish Balay #undef __FUNCT__ 20144a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 20153c7409f5SSatish Balay /*@C 2016f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2017d850072dSLois Curfman McInnes SNES options in the database. 20183c7409f5SSatish Balay 2019fee21e36SBarry Smith Collective on SNES 2020fee21e36SBarry Smith 2021c7afd0dbSLois Curfman McInnes Input Parameters: 2022c7afd0dbSLois Curfman McInnes + snes - the SNES context 2023c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2024c7afd0dbSLois Curfman McInnes 2025d850072dSLois Curfman McInnes Notes: 2026a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2027c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2028d850072dSLois Curfman McInnes 202936851e7fSLois Curfman McInnes Level: advanced 203036851e7fSLois Curfman McInnes 20313c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2032a86d99e1SLois Curfman McInnes 2033a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20343c7409f5SSatish Balay @*/ 2035dfbe8321SBarry Smith PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 20363c7409f5SSatish Balay { 2037dfbe8321SBarry Smith PetscErrorCode ierr; 20383c7409f5SSatish Balay 20393a40ed3dSBarry Smith PetscFunctionBegin; 20404482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2041639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 204294b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20433a40ed3dSBarry Smith PetscFunctionReturn(0); 20443c7409f5SSatish Balay } 20453c7409f5SSatish Balay 20464a2ae208SSatish Balay #undef __FUNCT__ 20474a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 20489ab63eb5SSatish Balay /*@C 20493c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20503c7409f5SSatish Balay SNES options in the database. 20513c7409f5SSatish Balay 2052c7afd0dbSLois Curfman McInnes Not Collective 2053c7afd0dbSLois Curfman McInnes 20543c7409f5SSatish Balay Input Parameter: 20553c7409f5SSatish Balay . snes - the SNES context 20563c7409f5SSatish Balay 20573c7409f5SSatish Balay Output Parameter: 20583c7409f5SSatish Balay . prefix - pointer to the prefix string used 20593c7409f5SSatish Balay 20609ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 20619ab63eb5SSatish Balay sufficient length to hold the prefix. 20629ab63eb5SSatish Balay 206336851e7fSLois Curfman McInnes Level: advanced 206436851e7fSLois Curfman McInnes 20653c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2066a86d99e1SLois Curfman McInnes 2067a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20683c7409f5SSatish Balay @*/ 2069dfbe8321SBarry Smith PetscErrorCode SNESGetOptionsPrefix(SNES snes,char *prefix[]) 20703c7409f5SSatish Balay { 2071dfbe8321SBarry Smith PetscErrorCode ierr; 20723c7409f5SSatish Balay 20733a40ed3dSBarry Smith PetscFunctionBegin; 20744482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2075639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 20763a40ed3dSBarry Smith PetscFunctionReturn(0); 20773c7409f5SSatish Balay } 20783c7409f5SSatish Balay 2079b2002411SLois Curfman McInnes 20804a2ae208SSatish Balay #undef __FUNCT__ 20814a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 20823cea93caSBarry Smith /*@C 20833cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 20843cea93caSBarry Smith 20857f6c08e0SMatthew Knepley Level: advanced 20863cea93caSBarry Smith @*/ 20876849ba73SBarry Smith PetscErrorCode SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2088b2002411SLois Curfman McInnes { 2089e2d1d2b7SBarry Smith char fullname[PETSC_MAX_PATH_LEN]; 2090dfbe8321SBarry Smith PetscErrorCode ierr; 2091b2002411SLois Curfman McInnes 2092b2002411SLois Curfman McInnes PetscFunctionBegin; 2093b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2094c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2095b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2096b2002411SLois Curfman McInnes } 2097da9b6338SBarry Smith 2098da9b6338SBarry Smith #undef __FUNCT__ 2099da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 2100dfbe8321SBarry Smith PetscErrorCode SNESTestLocalMin(SNES snes) 2101da9b6338SBarry Smith { 2102dfbe8321SBarry Smith PetscErrorCode ierr; 2103dfbe8321SBarry Smith int N,i,j; 2104da9b6338SBarry Smith Vec u,uh,fh; 2105da9b6338SBarry Smith PetscScalar value; 2106da9b6338SBarry Smith PetscReal norm; 2107da9b6338SBarry Smith 2108da9b6338SBarry Smith PetscFunctionBegin; 2109da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2110da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2111da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2112da9b6338SBarry Smith 2113da9b6338SBarry Smith /* currently only works for sequential */ 2114da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2115da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2116da9b6338SBarry Smith for (i=0; i<N; i++) { 2117da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 2118da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %d\n",i);CHKERRQ(ierr); 2119da9b6338SBarry Smith for (j=-10; j<11; j++) { 2120ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2121da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 21223ab0aad5SBarry Smith ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2123da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2124da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %d %18.16e\n",j,norm);CHKERRQ(ierr); 2125da9b6338SBarry Smith value = -value; 2126da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2127da9b6338SBarry Smith } 2128da9b6338SBarry Smith } 2129da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2130da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2131da9b6338SBarry Smith PetscFunctionReturn(0); 2132da9b6338SBarry Smith } 2133