173f4d377SMatthew Knepley /*$Id: snes.c,v 1.235 2001/08/21 21:03:49 bsmith Exp $*/ 29b94acceSBarry Smith 3e090d566SSatish Balay #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 49b94acceSBarry Smith 54c49b128SBarry Smith PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 68ba1e511SMatthew Knepley PetscFList SNESList = PETSC_NULL; 78ba1e511SMatthew Knepley 88ba1e511SMatthew Knepley /* Logging support */ 943c77886SBarry Smith int SNES_COOKIE = 0; 1043c77886SBarry Smith int SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0; 11a09944afSBarry Smith 12a09944afSBarry Smith #undef __FUNCT__ 134a2ae208SSatish Balay #define __FUNCT__ "SNESView" 147e2c5f70SBarry Smith /*@C 159b94acceSBarry Smith SNESView - Prints the SNES data structure. 169b94acceSBarry Smith 174c49b128SBarry Smith Collective on SNES 18fee21e36SBarry Smith 19c7afd0dbSLois Curfman McInnes Input Parameters: 20c7afd0dbSLois Curfman McInnes + SNES - the SNES context 21c7afd0dbSLois Curfman McInnes - viewer - visualization context 22c7afd0dbSLois Curfman McInnes 239b94acceSBarry Smith Options Database Key: 24c8a8ba5cSLois Curfman McInnes . -snes_view - Calls SNESView() at end of SNESSolve() 259b94acceSBarry Smith 269b94acceSBarry Smith Notes: 279b94acceSBarry Smith The available visualization contexts include 28b0a32e0cSBarry Smith + PETSC_VIEWER_STDOUT_SELF - standard output (default) 29b0a32e0cSBarry Smith - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 30c8a8ba5cSLois Curfman McInnes output where only the first processor opens 31c8a8ba5cSLois Curfman McInnes the file. All other processors send their 32c8a8ba5cSLois Curfman McInnes data to the first processor to print. 339b94acceSBarry Smith 343e081fefSLois Curfman McInnes The user can open an alternative visualization context with 35b0a32e0cSBarry Smith PetscViewerASCIIOpen() - output to a specified file. 369b94acceSBarry Smith 3736851e7fSLois Curfman McInnes Level: beginner 3836851e7fSLois Curfman McInnes 399b94acceSBarry Smith .keywords: SNES, view 409b94acceSBarry Smith 41b0a32e0cSBarry Smith .seealso: PetscViewerASCIIOpen() 429b94acceSBarry Smith @*/ 43b0a32e0cSBarry Smith int SNESView(SNES snes,PetscViewer viewer) 449b94acceSBarry Smith { 459b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 469b94acceSBarry Smith int ierr; 4794b7f48cSBarry Smith KSP ksp; 48454a90a3SBarry Smith char *type; 496831982aSBarry Smith PetscTruth isascii,isstring; 509b94acceSBarry Smith 513a40ed3dSBarry Smith PetscFunctionBegin; 524482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 53b0a32e0cSBarry Smith if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm); 544482741eSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 55c9780b6fSBarry Smith PetscCheckSameComm(snes,1,viewer,2); 5674679c65SBarry Smith 57b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 58b0a32e0cSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 590f5bd95cSBarry Smith if (isascii) { 603a7fca6bSBarry Smith if (snes->prefix) { 613a7fca6bSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr); 623a7fca6bSBarry Smith } else { 63b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 643a7fca6bSBarry Smith } 65454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 66454a90a3SBarry Smith if (type) { 67b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 68184914b5SBarry Smith } else { 69b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 70184914b5SBarry Smith } 710ef38995SBarry Smith if (snes->view) { 72b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 730ef38995SBarry Smith ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 74b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 750ef38995SBarry Smith } 76b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 7777d8c4bbSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n", 7877d8c4bbSBarry Smith snes->rtol,snes->atol,snes->xtol);CHKERRQ(ierr); 79b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr); 80b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr); 819b94acceSBarry Smith if (snes->ksp_ewconv) { 829b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 839b94acceSBarry Smith if (kctx) { 84b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr); 85b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 86b0a32e0cSBarry Smith ierr = PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 879b94acceSBarry Smith } 889b94acceSBarry Smith } 890f5bd95cSBarry Smith } else if (isstring) { 90454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 91b0a32e0cSBarry Smith ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 9219bcc07fSBarry Smith } 9394b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 94b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 9594b7f48cSBarry Smith ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 96b0a32e0cSBarry Smith ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 973a40ed3dSBarry Smith PetscFunctionReturn(0); 989b94acceSBarry Smith } 999b94acceSBarry Smith 10076b2cf59SMatthew Knepley /* 10176b2cf59SMatthew Knepley We retain a list of functions that also take SNES command 10276b2cf59SMatthew Knepley line options. These are called at the end SNESSetFromOptions() 10376b2cf59SMatthew Knepley */ 10476b2cf59SMatthew Knepley #define MAXSETFROMOPTIONS 5 10576b2cf59SMatthew Knepley static int numberofsetfromoptions; 10676b2cf59SMatthew Knepley static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 10776b2cf59SMatthew Knepley 108e74ef692SMatthew Knepley #undef __FUNCT__ 109e74ef692SMatthew Knepley #define __FUNCT__ "SNESAddOptionsChecker" 110ac226902SBarry Smith /*@C 11176b2cf59SMatthew Knepley SNESAddOptionsChecker - Adds an additional function to check for SNES options. 11276b2cf59SMatthew Knepley 11376b2cf59SMatthew Knepley Not Collective 11476b2cf59SMatthew Knepley 11576b2cf59SMatthew Knepley Input Parameter: 11676b2cf59SMatthew Knepley . snescheck - function that checks for options 11776b2cf59SMatthew Knepley 11876b2cf59SMatthew Knepley Level: developer 11976b2cf59SMatthew Knepley 12076b2cf59SMatthew Knepley .seealso: SNESSetFromOptions() 12176b2cf59SMatthew Knepley @*/ 12276b2cf59SMatthew Knepley int SNESAddOptionsChecker(int (*snescheck)(SNES)) 12376b2cf59SMatthew Knepley { 12476b2cf59SMatthew Knepley PetscFunctionBegin; 12576b2cf59SMatthew Knepley if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 12676b2cf59SMatthew Knepley SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %d allowed", MAXSETFROMOPTIONS); 12776b2cf59SMatthew Knepley } 12876b2cf59SMatthew Knepley 12976b2cf59SMatthew Knepley othersetfromoptions[numberofsetfromoptions++] = snescheck; 13076b2cf59SMatthew Knepley PetscFunctionReturn(0); 13176b2cf59SMatthew Knepley } 13276b2cf59SMatthew Knepley 1334a2ae208SSatish Balay #undef __FUNCT__ 1344a2ae208SSatish Balay #define __FUNCT__ "SNESSetFromOptions" 1359b94acceSBarry Smith /*@ 13694b7f48cSBarry Smith SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 1379b94acceSBarry Smith 138c7afd0dbSLois Curfman McInnes Collective on SNES 139c7afd0dbSLois Curfman McInnes 1409b94acceSBarry Smith Input Parameter: 1419b94acceSBarry Smith . snes - the SNES context 1429b94acceSBarry Smith 14336851e7fSLois Curfman McInnes Options Database Keys: 1446831982aSBarry Smith + -snes_type <type> - ls, tr, umls, umtr, test 14582738288SBarry Smith . -snes_stol - convergence tolerance in terms of the norm 14682738288SBarry Smith of the change in the solution between steps 147b39c3a46SLois Curfman McInnes . -snes_atol <atol> - absolute tolerance of residual norm 148b39c3a46SLois Curfman McInnes . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 149b39c3a46SLois Curfman McInnes . -snes_max_it <max_it> - maximum number of iterations 150b39c3a46SLois Curfman McInnes . -snes_max_funcs <max_funcs> - maximum number of function evaluations 15150ffb88aSMatthew Knepley . -snes_max_fail <max_fail> - maximum number of failures 152b39c3a46SLois Curfman McInnes . -snes_trtol <trtol> - trust region tolerance 15382738288SBarry Smith . -snes_no_convergence_test - skip convergence test in nonlinear or minimization 15482738288SBarry Smith solver; hence iterations will continue until max_it 1551fbbfb26SLois Curfman McInnes or some other criterion is reached. Saves expense 15682738288SBarry Smith of convergence test 15782738288SBarry Smith . -snes_monitor - prints residual norm at each iteration 158d132466eSBarry Smith . -snes_vecmonitor - plots solution at each iteration 159d132466eSBarry Smith . -snes_vecmonitor_update - plots update to solution at each iteration 16082738288SBarry Smith . -snes_xmonitor - plots residual norm at each iteration 161e24b481bSBarry Smith . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 1625968eb51SBarry Smith . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 1635968eb51SBarry Smith - -snes_print_converged_reason - print the reason for convergence/divergence after each solve 16482738288SBarry Smith 16582738288SBarry Smith Options Database for Eisenstat-Walker method: 1664b27c08aSLois Curfman McInnes + -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence 1674b27c08aSLois Curfman McInnes . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 16836851e7fSLois Curfman McInnes . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 16936851e7fSLois Curfman McInnes . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 17036851e7fSLois Curfman McInnes . -snes_ksp_ew_gamma <gamma> - Sets gamma 17136851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha <alpha> - Sets alpha 17236851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 17336851e7fSLois Curfman McInnes - -snes_ksp_ew_threshold <threshold> - Sets threshold 17482738288SBarry Smith 17511ca99fdSLois Curfman McInnes Notes: 17611ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 17711ca99fdSLois Curfman McInnes the users manual. 17883e2fdc7SBarry Smith 17936851e7fSLois Curfman McInnes Level: beginner 18036851e7fSLois Curfman McInnes 1819b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1829b94acceSBarry Smith 18369ed319cSSatish Balay .seealso: SNESSetOptionsPrefix() 1849b94acceSBarry Smith @*/ 1859b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1869b94acceSBarry Smith { 18794b7f48cSBarry Smith KSP ksp; 188186905e3SBarry Smith SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 189f1af5d2fSBarry Smith PetscTruth flg; 19076b2cf59SMatthew Knepley int ierr, 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 @*/ 2959b94acceSBarry Smith int 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 @*/ 3239b94acceSBarry Smith int 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 @*/ 3639b94acceSBarry Smith int 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 @*/ 39287828ca2SBarry Smith int 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 @*/ 4229b94acceSBarry Smith int 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 @*/ 44750ffb88aSMatthew Knepley int 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 @*/ 47350ffb88aSMatthew Knepley int 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 @*/ 50386bddb7dSBarry Smith int 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 @*/ 53694b7f48cSBarry Smith int 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" 547454a90a3SBarry Smith static int SNESPublish_Petsc(PetscObject obj) 548e24b481bSBarry Smith { 549aa482453SBarry Smith #if defined(PETSC_HAVE_AMS) 550454a90a3SBarry Smith SNES v = (SNES) obj; 551e24b481bSBarry Smith int ierr; 55243d6d2cbSBarry Smith #endif 553e24b481bSBarry Smith 554e24b481bSBarry Smith PetscFunctionBegin; 555e24b481bSBarry Smith 55643d6d2cbSBarry Smith #if defined(PETSC_HAVE_AMS) 557e24b481bSBarry Smith /* if it is already published then return */ 558e24b481bSBarry Smith if (v->amem >=0) PetscFunctionReturn(0); 559e24b481bSBarry Smith 560454a90a3SBarry Smith ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr); 561e24b481bSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ, 562e24b481bSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 563e24b481bSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ, 564e24b481bSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 565454a90a3SBarry Smith ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr); 566433994e6SBarry Smith #endif 567e24b481bSBarry Smith PetscFunctionReturn(0); 568e24b481bSBarry Smith } 569e24b481bSBarry Smith 5709b94acceSBarry Smith /* -----------------------------------------------------------*/ 5714a2ae208SSatish Balay #undef __FUNCT__ 5724a2ae208SSatish Balay #define __FUNCT__ "SNESCreate" 5739b94acceSBarry Smith /*@C 5749b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5759b94acceSBarry Smith 576c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 577c7afd0dbSLois Curfman McInnes 578c7afd0dbSLois Curfman McInnes Input Parameters: 579c7afd0dbSLois Curfman McInnes + comm - MPI communicator 5809b94acceSBarry Smith 5819b94acceSBarry Smith Output Parameter: 5829b94acceSBarry Smith . outsnes - the new SNES context 5839b94acceSBarry Smith 584c7afd0dbSLois Curfman McInnes Options Database Keys: 585c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 586c7afd0dbSLois Curfman McInnes and no preconditioning matrix 587c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 588c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 589c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 590c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 591c1f60f51SBarry Smith 59236851e7fSLois Curfman McInnes Level: beginner 59336851e7fSLois Curfman McInnes 5949b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5959b94acceSBarry Smith 5964b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES 5979b94acceSBarry Smith @*/ 5984b27c08aSLois Curfman McInnes int SNESCreate(MPI_Comm comm,SNES *outsnes) 5999b94acceSBarry Smith { 6009b94acceSBarry Smith int ierr; 6019b94acceSBarry Smith SNES snes; 6029b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 60337fcc0dbSBarry Smith 6043a40ed3dSBarry Smith PetscFunctionBegin; 6054482741eSBarry Smith PetscValidPointer(outsnes,1); 6068ba1e511SMatthew Knepley *outsnes = PETSC_NULL; 6078ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 6088ba1e511SMatthew Knepley ierr = SNESInitializePackage(PETSC_NULL); CHKERRQ(ierr); 6098ba1e511SMatthew Knepley #endif 6108ba1e511SMatthew Knepley 6113f1db9ecSBarry Smith PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView); 612b0a32e0cSBarry Smith PetscLogObjectCreate(snes); 613e24b481bSBarry Smith snes->bops->publish = SNESPublish_Petsc; 6149b94acceSBarry Smith snes->max_its = 50; 6159750a799SBarry Smith snes->max_funcs = 10000; 6169b94acceSBarry Smith snes->norm = 0.0; 617b4874afaSBarry Smith snes->rtol = 1.e-8; 618b4874afaSBarry Smith snes->ttol = 0.0; 619b4874afaSBarry Smith snes->atol = 1.e-50; 6209b94acceSBarry Smith snes->xtol = 1.e-8; 6214b27c08aSLois Curfman McInnes snes->deltatol = 1.e-12; 6229b94acceSBarry Smith snes->nfuncs = 0; 62350ffb88aSMatthew Knepley snes->numFailures = 0; 62450ffb88aSMatthew Knepley snes->maxFailures = 1; 6257a00f4a9SLois Curfman McInnes snes->linear_its = 0; 626639f9d9dSBarry Smith snes->numbermonitors = 0; 6279b94acceSBarry Smith snes->data = 0; 6289b94acceSBarry Smith snes->view = 0; 62982bf6240SBarry Smith snes->setupcalled = 0; 630186905e3SBarry Smith snes->ksp_ewconv = PETSC_FALSE; 6316f24a144SLois Curfman McInnes snes->vwork = 0; 6326f24a144SLois Curfman McInnes snes->nwork = 0; 633758f92a0SBarry Smith snes->conv_hist_len = 0; 634758f92a0SBarry Smith snes->conv_hist_max = 0; 635758f92a0SBarry Smith snes->conv_hist = PETSC_NULL; 636758f92a0SBarry Smith snes->conv_hist_its = PETSC_NULL; 637758f92a0SBarry Smith snes->conv_hist_reset = PETSC_TRUE; 638184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 6399b94acceSBarry Smith 6409b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 641b0a32e0cSBarry Smith ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 642b0a32e0cSBarry Smith PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 6439b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6449b94acceSBarry Smith kctx->version = 2; 6459b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6469b94acceSBarry Smith this was too large for some test cases */ 6479b94acceSBarry Smith kctx->rtol_last = 0; 6489b94acceSBarry Smith kctx->rtol_max = .9; 6499b94acceSBarry Smith kctx->gamma = 1.0; 6509b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6519b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6529b94acceSBarry Smith kctx->threshold = .1; 6539b94acceSBarry Smith kctx->lresid_last = 0; 6549b94acceSBarry Smith kctx->norm_last = 0; 6559b94acceSBarry Smith 65694b7f48cSBarry Smith ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr); 65794b7f48cSBarry Smith PetscLogObjectParent(snes,snes->ksp) 6585334005bSBarry Smith 6599b94acceSBarry Smith *outsnes = snes; 66000036973SBarry Smith ierr = PetscPublishAll(snes);CHKERRQ(ierr); 6613a40ed3dSBarry Smith PetscFunctionReturn(0); 6629b94acceSBarry Smith } 6639b94acceSBarry Smith 6644a2ae208SSatish Balay #undef __FUNCT__ 6654a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction" 6669b94acceSBarry Smith /*@C 6679b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6689b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6699b94acceSBarry Smith equations. 6709b94acceSBarry Smith 671fee21e36SBarry Smith Collective on SNES 672fee21e36SBarry Smith 673c7afd0dbSLois Curfman McInnes Input Parameters: 674c7afd0dbSLois Curfman McInnes + snes - the SNES context 675c7afd0dbSLois Curfman McInnes . func - function evaluation routine 676c7afd0dbSLois Curfman McInnes . r - vector to store function value 677c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 678c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6799b94acceSBarry Smith 680c7afd0dbSLois Curfman McInnes Calling sequence of func: 6818d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 682c7afd0dbSLois Curfman McInnes 683313e4042SLois Curfman McInnes . f - function vector 684c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 6859b94acceSBarry Smith 6869b94acceSBarry Smith Notes: 6879b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6889b94acceSBarry Smith $ f'(x) x = -f(x), 689c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 6909b94acceSBarry Smith 69136851e7fSLois Curfman McInnes Level: beginner 69236851e7fSLois Curfman McInnes 6939b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6949b94acceSBarry Smith 695a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6969b94acceSBarry Smith @*/ 6975334005bSBarry Smith int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 6989b94acceSBarry Smith { 6993a40ed3dSBarry Smith PetscFunctionBegin; 7004482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7014482741eSBarry Smith PetscValidHeaderSpecific(r,VEC_COOKIE,2); 702c9780b6fSBarry Smith PetscCheckSameComm(snes,1,r,2); 703184914b5SBarry Smith 7049b94acceSBarry Smith snes->computefunction = func; 7059b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7069b94acceSBarry Smith snes->funP = ctx; 7073a40ed3dSBarry Smith PetscFunctionReturn(0); 7089b94acceSBarry Smith } 7099b94acceSBarry Smith 710*3ab0aad5SBarry Smith /* --------------------------------------------------------------- */ 711*3ab0aad5SBarry Smith #undef __FUNCT__ 712*3ab0aad5SBarry Smith #define __FUNCT__ "SNESSetRhs" 713*3ab0aad5SBarry Smith /*@C 714*3ab0aad5SBarry Smith SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set 715*3ab0aad5SBarry Smith it assumes a zero right hand side. 716*3ab0aad5SBarry Smith 717*3ab0aad5SBarry Smith Collective on SNES 718*3ab0aad5SBarry Smith 719*3ab0aad5SBarry Smith Input Parameters: 720*3ab0aad5SBarry Smith + snes - the SNES context 721*3ab0aad5SBarry Smith - rhs - the right hand side vector or PETSC_NULL for a zero right hand side 722*3ab0aad5SBarry Smith 723*3ab0aad5SBarry Smith Level: intermediate 724*3ab0aad5SBarry Smith 725*3ab0aad5SBarry Smith .keywords: SNES, nonlinear, set, function, right hand side 726*3ab0aad5SBarry Smith 727*3ab0aad5SBarry Smith .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 728*3ab0aad5SBarry Smith @*/ 729*3ab0aad5SBarry Smith int SNESSetRhs(SNES snes,Vec rhs) 730*3ab0aad5SBarry Smith { 731*3ab0aad5SBarry Smith int ierr; 732*3ab0aad5SBarry Smith 733*3ab0aad5SBarry Smith PetscFunctionBegin; 734*3ab0aad5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 735*3ab0aad5SBarry Smith if (rhs) { 736*3ab0aad5SBarry Smith PetscValidHeaderSpecific(rhs,VEC_COOKIE,2); 737*3ab0aad5SBarry Smith PetscCheckSameComm(snes,1,rhs,2); 738*3ab0aad5SBarry Smith ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr); 739*3ab0aad5SBarry Smith } 740*3ab0aad5SBarry Smith if (snes->afine) { 741*3ab0aad5SBarry Smith ierr = VecDestroy(snes->afine);CHKERRQ(ierr); 742*3ab0aad5SBarry Smith } 743*3ab0aad5SBarry Smith snes->afine = rhs; 744*3ab0aad5SBarry Smith PetscFunctionReturn(0); 745*3ab0aad5SBarry Smith } 746*3ab0aad5SBarry Smith 7474a2ae208SSatish Balay #undef __FUNCT__ 7484a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction" 7499b94acceSBarry Smith /*@ 75036851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7519b94acceSBarry Smith SNESSetFunction(). 7529b94acceSBarry Smith 753c7afd0dbSLois Curfman McInnes Collective on SNES 754c7afd0dbSLois Curfman McInnes 7559b94acceSBarry Smith Input Parameters: 756c7afd0dbSLois Curfman McInnes + snes - the SNES context 757c7afd0dbSLois Curfman McInnes - x - input vector 7589b94acceSBarry Smith 7599b94acceSBarry Smith Output Parameter: 7603638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7619b94acceSBarry Smith 7621bffabb2SLois Curfman McInnes Notes: 76336851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 76436851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 76536851e7fSLois Curfman McInnes themselves. 76636851e7fSLois Curfman McInnes 76736851e7fSLois Curfman McInnes Level: developer 76836851e7fSLois Curfman McInnes 7699b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 7709b94acceSBarry Smith 771a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 7729b94acceSBarry Smith @*/ 7739b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x,Vec y) 7749b94acceSBarry Smith { 7759b94acceSBarry Smith int ierr; 7769b94acceSBarry Smith 7773a40ed3dSBarry Smith PetscFunctionBegin; 7784482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7794482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 7804482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,3); 781c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 782c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,3); 783184914b5SBarry Smith 784d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 785d64ed03dSBarry Smith PetscStackPush("SNES user function"); 7869b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 787d64ed03dSBarry Smith PetscStackPop; 788*3ab0aad5SBarry Smith if (snes->afine) { 789*3ab0aad5SBarry Smith PetscScalar mone = -1.0; 790*3ab0aad5SBarry Smith ierr = VecAXPY(&mone,snes->afine,y);CHKERRQ(ierr); 791*3ab0aad5SBarry Smith } 792ae3c334cSLois Curfman McInnes snes->nfuncs++; 793d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 7943a40ed3dSBarry Smith PetscFunctionReturn(0); 7959b94acceSBarry Smith } 7969b94acceSBarry Smith 7974a2ae208SSatish Balay #undef __FUNCT__ 7984a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian" 79962fef451SLois Curfman McInnes /*@ 80062fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 80162fef451SLois Curfman McInnes set with SNESSetJacobian(). 80262fef451SLois Curfman McInnes 803c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 804c7afd0dbSLois Curfman McInnes 80562fef451SLois Curfman McInnes Input Parameters: 806c7afd0dbSLois Curfman McInnes + snes - the SNES context 807c7afd0dbSLois Curfman McInnes - x - input vector 80862fef451SLois Curfman McInnes 80962fef451SLois Curfman McInnes Output Parameters: 810c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 81162fef451SLois Curfman McInnes . B - optional preconditioning matrix 812c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 813fee21e36SBarry Smith 81462fef451SLois Curfman McInnes Notes: 81562fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 81662fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 81762fef451SLois Curfman McInnes 81894b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 819dc5a77f8SLois Curfman McInnes flag parameter. 82062fef451SLois Curfman McInnes 82136851e7fSLois Curfman McInnes Level: developer 82236851e7fSLois Curfman McInnes 82362fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 82462fef451SLois Curfman McInnes 82594b7f48cSBarry Smith .seealso: SNESSetJacobian(), KSPSetOperators() 82662fef451SLois Curfman McInnes @*/ 8279b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8289b94acceSBarry Smith { 8299b94acceSBarry Smith int ierr; 8303a40ed3dSBarry Smith 8313a40ed3dSBarry Smith PetscFunctionBegin; 8324482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8334482741eSBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE,2); 8344482741eSBarry Smith PetscValidPointer(flg,5); 835c9780b6fSBarry Smith PetscCheckSameComm(snes,1,X,2); 8363a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 837d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 838c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 839d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 8409b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 841d64ed03dSBarry Smith PetscStackPop; 842d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 8436d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 8444482741eSBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 8454482741eSBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 8463a40ed3dSBarry Smith PetscFunctionReturn(0); 8479b94acceSBarry Smith } 8489b94acceSBarry Smith 8494a2ae208SSatish Balay #undef __FUNCT__ 8504a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian" 8519b94acceSBarry Smith /*@C 8529b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 853044dda88SLois Curfman McInnes location to store the matrix. 8549b94acceSBarry Smith 855c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 856c7afd0dbSLois Curfman McInnes 8579b94acceSBarry Smith Input Parameters: 858c7afd0dbSLois Curfman McInnes + snes - the SNES context 8599b94acceSBarry Smith . A - Jacobian matrix 8609b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 8619b94acceSBarry Smith . func - Jacobian evaluation routine 862c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 8632cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 8649b94acceSBarry Smith 8659b94acceSBarry Smith Calling sequence of func: 8668d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 8679b94acceSBarry Smith 868c7afd0dbSLois Curfman McInnes + x - input vector 8699b94acceSBarry Smith . A - Jacobian matrix 8709b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 871ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 87294b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 873c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 8749b94acceSBarry Smith 8759b94acceSBarry Smith Notes: 87694b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 8772cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 878ac21db08SLois Curfman McInnes 879ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 8809b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 8819b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 8829b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 8839b94acceSBarry Smith throughout the global iterations. 8849b94acceSBarry Smith 88536851e7fSLois Curfman McInnes Level: beginner 88636851e7fSLois Curfman McInnes 8879b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 8889b94acceSBarry Smith 88994b7f48cSBarry Smith .seealso: KSPSetOperators(), SNESSetFunction() 8909b94acceSBarry Smith @*/ 891454a90a3SBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 8929b94acceSBarry Smith { 8933a7fca6bSBarry Smith int ierr; 8943a7fca6bSBarry Smith 8953a40ed3dSBarry Smith PetscFunctionBegin; 8964482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8974482741eSBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 8984482741eSBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 899c9780b6fSBarry Smith if (A) PetscCheckSameComm(snes,1,A,2); 900c9780b6fSBarry Smith if (B) PetscCheckSameComm(snes,1,B,2); 9013a7fca6bSBarry Smith if (func) snes->computejacobian = func; 9023a7fca6bSBarry Smith if (ctx) snes->jacP = ctx; 9033a7fca6bSBarry Smith if (A) { 9043a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 9059b94acceSBarry Smith snes->jacobian = A; 9063a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 9073a7fca6bSBarry Smith } 9083a7fca6bSBarry Smith if (B) { 9093a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 9109b94acceSBarry Smith snes->jacobian_pre = B; 9113a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 9123a7fca6bSBarry Smith } 9133a40ed3dSBarry Smith PetscFunctionReturn(0); 9149b94acceSBarry Smith } 91562fef451SLois Curfman McInnes 9164a2ae208SSatish Balay #undef __FUNCT__ 9174a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian" 918c2aafc4cSSatish Balay /*@C 919b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 920b4fd4287SBarry Smith provided context for evaluating the Jacobian. 921b4fd4287SBarry Smith 922c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 923c7afd0dbSLois Curfman McInnes 924b4fd4287SBarry Smith Input Parameter: 925b4fd4287SBarry Smith . snes - the nonlinear solver context 926b4fd4287SBarry Smith 927b4fd4287SBarry Smith Output Parameters: 928c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 929b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 93000036973SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 93100036973SBarry Smith - func - location to put Jacobian function (or PETSC_NULL) 932fee21e36SBarry Smith 93336851e7fSLois Curfman McInnes Level: advanced 93436851e7fSLois Curfman McInnes 935b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 936b4fd4287SBarry Smith @*/ 93700036973SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 938b4fd4287SBarry Smith { 9393a40ed3dSBarry Smith PetscFunctionBegin; 9404482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 941b4fd4287SBarry Smith if (A) *A = snes->jacobian; 942b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 943b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 94400036973SBarry Smith if (func) *func = snes->computejacobian; 9453a40ed3dSBarry Smith PetscFunctionReturn(0); 946b4fd4287SBarry Smith } 947b4fd4287SBarry Smith 9489b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 94945fc7adcSBarry Smith extern int SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 9509b94acceSBarry Smith 9514a2ae208SSatish Balay #undef __FUNCT__ 9524a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp" 9539b94acceSBarry Smith /*@ 9549b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 955272ac6f2SLois Curfman McInnes of a nonlinear solver. 9569b94acceSBarry Smith 957fee21e36SBarry Smith Collective on SNES 958fee21e36SBarry Smith 959c7afd0dbSLois Curfman McInnes Input Parameters: 960c7afd0dbSLois Curfman McInnes + snes - the SNES context 961c7afd0dbSLois Curfman McInnes - x - the solution vector 962c7afd0dbSLois Curfman McInnes 963272ac6f2SLois Curfman McInnes Notes: 964272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 965272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 966272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 967272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 968272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 969272ac6f2SLois Curfman McInnes 97036851e7fSLois Curfman McInnes Level: advanced 97136851e7fSLois Curfman McInnes 9729b94acceSBarry Smith .keywords: SNES, nonlinear, setup 9739b94acceSBarry Smith 9749b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 9759b94acceSBarry Smith @*/ 9768ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 9779b94acceSBarry Smith { 978f1af5d2fSBarry Smith int ierr; 9794b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 9803a40ed3dSBarry Smith 9813a40ed3dSBarry Smith PetscFunctionBegin; 9824482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 9834482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 984c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 9858ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 9869b94acceSBarry Smith 987b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 988c1f60f51SBarry Smith /* 989c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 990dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 991c1f60f51SBarry Smith */ 992c1f60f51SBarry Smith if (flg) { 993c1f60f51SBarry Smith Mat J; 9945a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 9955a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 9963a7fca6bSBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n"); 9973a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 9983a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 999c1f60f51SBarry Smith } 100045fc7adcSBarry Smith 1001b4014bb2SMatthew Knepley #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) 100245fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 100345fc7adcSBarry Smith if (flg) { 100445fc7adcSBarry Smith Mat J; 100545fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 100645fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 100745fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 100845fc7adcSBarry Smith } 100932a4b47aSMatthew Knepley #endif 101045fc7adcSBarry Smith 1011b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1012c1f60f51SBarry Smith /* 1013dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1014c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1015c1f60f51SBarry Smith */ 101631615d04SLois Curfman McInnes if (flg) { 1017272ac6f2SLois Curfman McInnes Mat J; 1018b5d62d44SBarry Smith KSP ksp; 101994b7f48cSBarry Smith PC pc; 1020f3ef73ceSBarry Smith 10215a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 10223a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 102393b98531SBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"); 10243a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 10253a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 10263a7fca6bSBarry Smith 1027f3ef73ceSBarry Smith /* force no preconditioner */ 102894b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1029b5d62d44SBarry Smith ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1030a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1031a9815358SBarry Smith if (!flg) { 1032f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1033272ac6f2SLois Curfman McInnes } 1034a9815358SBarry Smith } 1035f3ef73ceSBarry Smith 103629bbc08cSBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 103729bbc08cSBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 103829bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1039a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 104029bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1041a8c6a408SBarry Smith } 1042a703fe33SLois Curfman McInnes 1043a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 10444b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 10456831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 104694b7f48cSBarry Smith KSP ksp; 104794b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1048186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1049a703fe33SLois Curfman McInnes } 10504b27c08aSLois Curfman McInnes 1051a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 105282bf6240SBarry Smith snes->setupcalled = 1; 10533a40ed3dSBarry Smith PetscFunctionReturn(0); 10549b94acceSBarry Smith } 10559b94acceSBarry Smith 10564a2ae208SSatish Balay #undef __FUNCT__ 10574a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 10589b94acceSBarry Smith /*@C 10599b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 10609b94acceSBarry Smith with SNESCreate(). 10619b94acceSBarry Smith 1062c7afd0dbSLois Curfman McInnes Collective on SNES 1063c7afd0dbSLois Curfman McInnes 10649b94acceSBarry Smith Input Parameter: 10659b94acceSBarry Smith . snes - the SNES context 10669b94acceSBarry Smith 106736851e7fSLois Curfman McInnes Level: beginner 106836851e7fSLois Curfman McInnes 10699b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 10709b94acceSBarry Smith 107163a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 10729b94acceSBarry Smith @*/ 10739b94acceSBarry Smith int SNESDestroy(SNES snes) 10749b94acceSBarry Smith { 1075b8a78c4aSBarry Smith int i,ierr; 10763a40ed3dSBarry Smith 10773a40ed3dSBarry Smith PetscFunctionBegin; 10784482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 10793a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1080d4bb536fSBarry Smith 1081be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 10820f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1083be0abb6dSBarry Smith 1084e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1085606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 10863a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 10873a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1088*3ab0aad5SBarry Smith if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 108994b7f48cSBarry Smith ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1090522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1091b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1092b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1093b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1094b8a78c4aSBarry Smith } 1095b8a78c4aSBarry Smith } 1096b0a32e0cSBarry Smith PetscLogObjectDestroy((PetscObject)snes); 10970452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 10983a40ed3dSBarry Smith PetscFunctionReturn(0); 10999b94acceSBarry Smith } 11009b94acceSBarry Smith 11019b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 11029b94acceSBarry Smith 11034a2ae208SSatish Balay #undef __FUNCT__ 11044a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 11059b94acceSBarry Smith /*@ 1106d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 11079b94acceSBarry Smith 1108c7afd0dbSLois Curfman McInnes Collective on SNES 1109c7afd0dbSLois Curfman McInnes 11109b94acceSBarry Smith Input Parameters: 1111c7afd0dbSLois Curfman McInnes + snes - the SNES context 111233174efeSLois Curfman McInnes . atol - absolute convergence tolerance 111333174efeSLois Curfman McInnes . rtol - relative convergence tolerance 111433174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 111533174efeSLois Curfman McInnes of the change in the solution between steps 111633174efeSLois Curfman McInnes . maxit - maximum number of iterations 1117c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1118fee21e36SBarry Smith 111933174efeSLois Curfman McInnes Options Database Keys: 1120c7afd0dbSLois Curfman McInnes + -snes_atol <atol> - Sets atol 1121c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1122c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1123c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1124c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 11259b94acceSBarry Smith 1126d7a720efSLois Curfman McInnes Notes: 11279b94acceSBarry Smith The default maximum number of iterations is 50. 11289b94acceSBarry Smith The default maximum number of function evaluations is 1000. 11299b94acceSBarry Smith 113036851e7fSLois Curfman McInnes Level: intermediate 113136851e7fSLois Curfman McInnes 113233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 11339b94acceSBarry Smith 1134d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 11359b94acceSBarry Smith @*/ 1136329f5518SBarry Smith int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf) 11379b94acceSBarry Smith { 11383a40ed3dSBarry Smith PetscFunctionBegin; 11394482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1140d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1141d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1142d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1143d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1144d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11453a40ed3dSBarry Smith PetscFunctionReturn(0); 11469b94acceSBarry Smith } 11479b94acceSBarry Smith 11484a2ae208SSatish Balay #undef __FUNCT__ 11494a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11509b94acceSBarry Smith /*@ 115133174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 115233174efeSLois Curfman McInnes 1153c7afd0dbSLois Curfman McInnes Not Collective 1154c7afd0dbSLois Curfman McInnes 115533174efeSLois Curfman McInnes Input Parameters: 1156c7afd0dbSLois Curfman McInnes + snes - the SNES context 115733174efeSLois Curfman McInnes . atol - absolute convergence tolerance 115833174efeSLois Curfman McInnes . rtol - relative convergence tolerance 115933174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 116033174efeSLois Curfman McInnes of the change in the solution between steps 116133174efeSLois Curfman McInnes . maxit - maximum number of iterations 1162c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1163fee21e36SBarry Smith 116433174efeSLois Curfman McInnes Notes: 116533174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 116633174efeSLois Curfman McInnes 116736851e7fSLois Curfman McInnes Level: intermediate 116836851e7fSLois Curfman McInnes 116933174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 117033174efeSLois Curfman McInnes 117133174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 117233174efeSLois Curfman McInnes @*/ 1173329f5518SBarry Smith int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf) 117433174efeSLois Curfman McInnes { 11753a40ed3dSBarry Smith PetscFunctionBegin; 11764482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 117733174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 117833174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 117933174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 118033174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 118133174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 11823a40ed3dSBarry Smith PetscFunctionReturn(0); 118333174efeSLois Curfman McInnes } 118433174efeSLois Curfman McInnes 11854a2ae208SSatish Balay #undef __FUNCT__ 11864a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 118733174efeSLois Curfman McInnes /*@ 11889b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 11899b94acceSBarry Smith 1190fee21e36SBarry Smith Collective on SNES 1191fee21e36SBarry Smith 1192c7afd0dbSLois Curfman McInnes Input Parameters: 1193c7afd0dbSLois Curfman McInnes + snes - the SNES context 1194c7afd0dbSLois Curfman McInnes - tol - tolerance 1195c7afd0dbSLois Curfman McInnes 11969b94acceSBarry Smith Options Database Key: 1197c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 11989b94acceSBarry Smith 119936851e7fSLois Curfman McInnes Level: intermediate 120036851e7fSLois Curfman McInnes 12019b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 12029b94acceSBarry Smith 1203d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 12049b94acceSBarry Smith @*/ 1205329f5518SBarry Smith int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 12069b94acceSBarry Smith { 12073a40ed3dSBarry Smith PetscFunctionBegin; 12084482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 12099b94acceSBarry Smith snes->deltatol = tol; 12103a40ed3dSBarry Smith PetscFunctionReturn(0); 12119b94acceSBarry Smith } 12129b94acceSBarry Smith 1213df9fa365SBarry Smith /* 1214df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1215df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1216df9fa365SBarry Smith macros instead of functions 1217df9fa365SBarry Smith */ 12184a2ae208SSatish Balay #undef __FUNCT__ 12194a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 1220329f5518SBarry Smith int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx) 1221ce1608b8SBarry Smith { 1222ce1608b8SBarry Smith int ierr; 1223ce1608b8SBarry Smith 1224ce1608b8SBarry Smith PetscFunctionBegin; 12254482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1226ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1227ce1608b8SBarry Smith PetscFunctionReturn(0); 1228ce1608b8SBarry Smith } 1229ce1608b8SBarry Smith 12304a2ae208SSatish Balay #undef __FUNCT__ 12314a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 12321836bdbcSSatish Balay int SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1233df9fa365SBarry Smith { 1234df9fa365SBarry Smith int ierr; 1235df9fa365SBarry Smith 1236df9fa365SBarry Smith PetscFunctionBegin; 1237df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1238df9fa365SBarry Smith PetscFunctionReturn(0); 1239df9fa365SBarry Smith } 1240df9fa365SBarry Smith 12414a2ae208SSatish Balay #undef __FUNCT__ 12424a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 1243b0a32e0cSBarry Smith int SNESLGMonitorDestroy(PetscDrawLG draw) 1244df9fa365SBarry Smith { 1245df9fa365SBarry Smith int ierr; 1246df9fa365SBarry Smith 1247df9fa365SBarry Smith PetscFunctionBegin; 1248df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1249df9fa365SBarry Smith PetscFunctionReturn(0); 1250df9fa365SBarry Smith } 1251df9fa365SBarry Smith 12529b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12539b94acceSBarry Smith 12544a2ae208SSatish Balay #undef __FUNCT__ 12554a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12569b94acceSBarry Smith /*@C 12575cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12589b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12599b94acceSBarry Smith progress. 12609b94acceSBarry Smith 1261fee21e36SBarry Smith Collective on SNES 1262fee21e36SBarry Smith 1263c7afd0dbSLois Curfman McInnes Input Parameters: 1264c7afd0dbSLois Curfman McInnes + snes - the SNES context 1265c7afd0dbSLois Curfman McInnes . func - monitoring routine 1266b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1267b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desitre) 1268b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1269b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 12709b94acceSBarry Smith 1271c7afd0dbSLois Curfman McInnes Calling sequence of func: 1272329f5518SBarry Smith $ int func(SNES snes,int its, PetscReal norm,void *mctx) 1273c7afd0dbSLois Curfman McInnes 1274c7afd0dbSLois Curfman McInnes + snes - the SNES context 1275c7afd0dbSLois Curfman McInnes . its - iteration number 1276c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 127740a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 12789b94acceSBarry Smith 12799665c990SLois Curfman McInnes Options Database Keys: 1280c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1281c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1282c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1283c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1284c7afd0dbSLois Curfman McInnes been hardwired into a code by 1285c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1286c7afd0dbSLois Curfman McInnes does not cancel those set via 1287c7afd0dbSLois Curfman McInnes the options database. 12889665c990SLois Curfman McInnes 1289639f9d9dSBarry Smith Notes: 12906bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 12916bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 12926bc08f3fSLois Curfman McInnes order in which they were set. 1293639f9d9dSBarry Smith 129436851e7fSLois Curfman McInnes Level: intermediate 129536851e7fSLois Curfman McInnes 12969b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 12979b94acceSBarry Smith 12985cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 12999b94acceSBarry Smith @*/ 1300329f5518SBarry Smith int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *)) 13019b94acceSBarry Smith { 13023a40ed3dSBarry Smith PetscFunctionBegin; 13034482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1304639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 130529bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1306639f9d9dSBarry Smith } 1307639f9d9dSBarry Smith 1308639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1309b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1310639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 13113a40ed3dSBarry Smith PetscFunctionReturn(0); 13129b94acceSBarry Smith } 13139b94acceSBarry Smith 13144a2ae208SSatish Balay #undef __FUNCT__ 13154a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 13165cd90555SBarry Smith /*@C 13175cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 13185cd90555SBarry Smith 1319c7afd0dbSLois Curfman McInnes Collective on SNES 1320c7afd0dbSLois Curfman McInnes 13215cd90555SBarry Smith Input Parameters: 13225cd90555SBarry Smith . snes - the SNES context 13235cd90555SBarry Smith 13241a480d89SAdministrator Options Database Key: 1325c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1326c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1327c7afd0dbSLois Curfman McInnes set via the options database 13285cd90555SBarry Smith 13295cd90555SBarry Smith Notes: 13305cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 13315cd90555SBarry Smith 133236851e7fSLois Curfman McInnes Level: intermediate 133336851e7fSLois Curfman McInnes 13345cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 13355cd90555SBarry Smith 13365cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 13375cd90555SBarry Smith @*/ 13385cd90555SBarry Smith int SNESClearMonitor(SNES snes) 13395cd90555SBarry Smith { 13405cd90555SBarry Smith PetscFunctionBegin; 13414482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13425cd90555SBarry Smith snes->numbermonitors = 0; 13435cd90555SBarry Smith PetscFunctionReturn(0); 13445cd90555SBarry Smith } 13455cd90555SBarry Smith 13464a2ae208SSatish Balay #undef __FUNCT__ 13474a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13489b94acceSBarry Smith /*@C 13499b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13509b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13519b94acceSBarry Smith 1352fee21e36SBarry Smith Collective on SNES 1353fee21e36SBarry Smith 1354c7afd0dbSLois Curfman McInnes Input Parameters: 1355c7afd0dbSLois Curfman McInnes + snes - the SNES context 1356c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1357c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1358c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 13599b94acceSBarry Smith 1360c7afd0dbSLois Curfman McInnes Calling sequence of func: 1361329f5518SBarry Smith $ int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1362c7afd0dbSLois Curfman McInnes 1363c7afd0dbSLois Curfman McInnes + snes - the SNES context 1364c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1365184914b5SBarry Smith . reason - reason for convergence/divergence 1366c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 13674b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 13684b27c08aSLois Curfman McInnes - f - 2-norm of function 13699b94acceSBarry Smith 137036851e7fSLois Curfman McInnes Level: advanced 137136851e7fSLois Curfman McInnes 13729b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13739b94acceSBarry Smith 13744b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 13759b94acceSBarry Smith @*/ 1376329f5518SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 13779b94acceSBarry Smith { 13783a40ed3dSBarry Smith PetscFunctionBegin; 13794482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13809b94acceSBarry Smith (snes)->converged = func; 13819b94acceSBarry Smith (snes)->cnvP = cctx; 13823a40ed3dSBarry Smith PetscFunctionReturn(0); 13839b94acceSBarry Smith } 13849b94acceSBarry Smith 13854a2ae208SSatish Balay #undef __FUNCT__ 13864a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 1387184914b5SBarry Smith /*@C 1388184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1389184914b5SBarry Smith 1390184914b5SBarry Smith Not Collective 1391184914b5SBarry Smith 1392184914b5SBarry Smith Input Parameter: 1393184914b5SBarry Smith . snes - the SNES context 1394184914b5SBarry Smith 1395184914b5SBarry Smith Output Parameter: 1396e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1397184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1398184914b5SBarry Smith 1399184914b5SBarry Smith Level: intermediate 1400184914b5SBarry Smith 1401184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1402184914b5SBarry Smith 1403184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1404184914b5SBarry Smith 14054b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1406184914b5SBarry Smith @*/ 1407184914b5SBarry Smith int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1408184914b5SBarry Smith { 1409184914b5SBarry Smith PetscFunctionBegin; 14104482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14114482741eSBarry Smith PetscValidPointer(reason,2); 1412184914b5SBarry Smith *reason = snes->reason; 1413184914b5SBarry Smith PetscFunctionReturn(0); 1414184914b5SBarry Smith } 1415184914b5SBarry Smith 14164a2ae208SSatish Balay #undef __FUNCT__ 14174a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1418c9005455SLois Curfman McInnes /*@ 1419c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1420c9005455SLois Curfman McInnes 1421fee21e36SBarry Smith Collective on SNES 1422fee21e36SBarry Smith 1423c7afd0dbSLois Curfman McInnes Input Parameters: 1424c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1425c7afd0dbSLois Curfman McInnes . a - array to hold history 1426cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1427758f92a0SBarry Smith . na - size of a and its 142864731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1429758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1430c7afd0dbSLois Curfman McInnes 1431c9005455SLois Curfman McInnes Notes: 14324b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1433c9005455SLois Curfman McInnes at each step. 1434c9005455SLois Curfman McInnes 1435c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1436c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1437c9005455SLois Curfman McInnes during the section of code that is being timed. 1438c9005455SLois Curfman McInnes 143936851e7fSLois Curfman McInnes Level: intermediate 144036851e7fSLois Curfman McInnes 1441c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1442758f92a0SBarry Smith 144308405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1444758f92a0SBarry Smith 1445c9005455SLois Curfman McInnes @*/ 14461836bdbcSSatish Balay int SNESSetConvergenceHistory(SNES snes,PetscReal a[],int *its,int na,PetscTruth reset) 1447c9005455SLois Curfman McInnes { 14483a40ed3dSBarry Smith PetscFunctionBegin; 14494482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14504482741eSBarry Smith if (na) PetscValidScalarPointer(a,2); 1451c9005455SLois Curfman McInnes snes->conv_hist = a; 1452758f92a0SBarry Smith snes->conv_hist_its = its; 1453758f92a0SBarry Smith snes->conv_hist_max = na; 1454758f92a0SBarry Smith snes->conv_hist_reset = reset; 1455758f92a0SBarry Smith PetscFunctionReturn(0); 1456758f92a0SBarry Smith } 1457758f92a0SBarry Smith 14584a2ae208SSatish Balay #undef __FUNCT__ 14594a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 14600c4c9dddSBarry Smith /*@C 1461758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1462758f92a0SBarry Smith 1463758f92a0SBarry Smith Collective on SNES 1464758f92a0SBarry Smith 1465758f92a0SBarry Smith Input Parameter: 1466758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1467758f92a0SBarry Smith 1468758f92a0SBarry Smith Output Parameters: 1469758f92a0SBarry Smith . a - array to hold history 1470758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1471758f92a0SBarry Smith negative if not converged) for each solve. 1472758f92a0SBarry Smith - na - size of a and its 1473758f92a0SBarry Smith 1474758f92a0SBarry Smith Notes: 1475758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1476758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1477758f92a0SBarry Smith 1478758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1479758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1480758f92a0SBarry Smith during the section of code that is being timed. 1481758f92a0SBarry Smith 1482758f92a0SBarry Smith Level: intermediate 1483758f92a0SBarry Smith 1484758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1485758f92a0SBarry Smith 1486758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1487758f92a0SBarry Smith 1488758f92a0SBarry Smith @*/ 14891836bdbcSSatish Balay int SNESGetConvergenceHistory(SNES snes,PetscReal *a[],int *its[],int *na) 1490758f92a0SBarry Smith { 1491758f92a0SBarry Smith PetscFunctionBegin; 14924482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1493758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1494758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1495758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 14963a40ed3dSBarry Smith PetscFunctionReturn(0); 1497c9005455SLois Curfman McInnes } 1498c9005455SLois Curfman McInnes 1499e74ef692SMatthew Knepley #undef __FUNCT__ 1500e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC" 1501ac226902SBarry Smith /*@C 150276b2cf59SMatthew Knepley SNESSetRhsBC - Sets the function which applies boundary conditions 150376b2cf59SMatthew Knepley to the Rhs of each system. 150476b2cf59SMatthew Knepley 150576b2cf59SMatthew Knepley Collective on SNES 150676b2cf59SMatthew Knepley 150776b2cf59SMatthew Knepley Input Parameters: 150876b2cf59SMatthew Knepley . snes - The nonlinear solver context 150976b2cf59SMatthew Knepley . func - The function 151076b2cf59SMatthew Knepley 151176b2cf59SMatthew Knepley Calling sequence of func: 151276b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx); 151376b2cf59SMatthew Knepley 151476b2cf59SMatthew Knepley . rhs - The current rhs vector 151576b2cf59SMatthew Knepley . ctx - The user-context 151676b2cf59SMatthew Knepley 151776b2cf59SMatthew Knepley Level: intermediate 151876b2cf59SMatthew Knepley 151976b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 152076b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate() 152176b2cf59SMatthew Knepley @*/ 152276b2cf59SMatthew Knepley int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *)) 152376b2cf59SMatthew Knepley { 152476b2cf59SMatthew Knepley PetscFunctionBegin; 15254482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 152676b2cf59SMatthew Knepley snes->applyrhsbc = func; 152776b2cf59SMatthew Knepley PetscFunctionReturn(0); 152876b2cf59SMatthew Knepley } 152976b2cf59SMatthew Knepley 1530e74ef692SMatthew Knepley #undef __FUNCT__ 1531e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC" 153276b2cf59SMatthew Knepley /*@ 153376b2cf59SMatthew Knepley SNESDefaultRhsBC - The default boundary condition function which does nothing. 153476b2cf59SMatthew Knepley 153576b2cf59SMatthew Knepley Not collective 153676b2cf59SMatthew Knepley 153776b2cf59SMatthew Knepley Input Parameters: 153876b2cf59SMatthew Knepley . snes - The nonlinear solver context 153976b2cf59SMatthew Knepley . rhs - The Rhs 154076b2cf59SMatthew Knepley . ctx - The user-context 154176b2cf59SMatthew Knepley 154276b2cf59SMatthew Knepley Level: beginner 154376b2cf59SMatthew Knepley 154476b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 154576b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate() 154676b2cf59SMatthew Knepley @*/ 154776b2cf59SMatthew Knepley int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx) 154876b2cf59SMatthew Knepley { 154976b2cf59SMatthew Knepley PetscFunctionBegin; 155076b2cf59SMatthew Knepley PetscFunctionReturn(0); 155176b2cf59SMatthew Knepley } 155276b2cf59SMatthew Knepley 1553e74ef692SMatthew Knepley #undef __FUNCT__ 1554e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC" 1555ac226902SBarry Smith /*@C 155676b2cf59SMatthew Knepley SNESSetSolutionBC - Sets the function which applies boundary conditions 155776b2cf59SMatthew Knepley to the solution of each system. 155876b2cf59SMatthew Knepley 155976b2cf59SMatthew Knepley Collective on SNES 156076b2cf59SMatthew Knepley 156176b2cf59SMatthew Knepley Input Parameters: 156276b2cf59SMatthew Knepley . snes - The nonlinear solver context 156376b2cf59SMatthew Knepley . func - The function 156476b2cf59SMatthew Knepley 156576b2cf59SMatthew Knepley Calling sequence of func: 15669d06fb6bSMatthew Knepley . func (SNES snes, Vec rsol, void *ctx); 156776b2cf59SMatthew Knepley 156876b2cf59SMatthew Knepley . sol - The current solution vector 156976b2cf59SMatthew Knepley . ctx - The user-context 157076b2cf59SMatthew Knepley 157176b2cf59SMatthew Knepley Level: intermediate 157276b2cf59SMatthew Knepley 157376b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 157476b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate() 157576b2cf59SMatthew Knepley @*/ 157676b2cf59SMatthew Knepley int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *)) 157776b2cf59SMatthew Knepley { 157876b2cf59SMatthew Knepley PetscFunctionBegin; 15794482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 158076b2cf59SMatthew Knepley snes->applysolbc = func; 158176b2cf59SMatthew Knepley PetscFunctionReturn(0); 158276b2cf59SMatthew Knepley } 158376b2cf59SMatthew Knepley 1584e74ef692SMatthew Knepley #undef __FUNCT__ 1585e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC" 158676b2cf59SMatthew Knepley /*@ 158776b2cf59SMatthew Knepley SNESDefaultSolutionBC - The default boundary condition function which does nothing. 158876b2cf59SMatthew Knepley 158976b2cf59SMatthew Knepley Not collective 159076b2cf59SMatthew Knepley 159176b2cf59SMatthew Knepley Input Parameters: 159276b2cf59SMatthew Knepley . snes - The nonlinear solver context 159376b2cf59SMatthew Knepley . sol - The solution 159476b2cf59SMatthew Knepley . ctx - The user-context 159576b2cf59SMatthew Knepley 159676b2cf59SMatthew Knepley Level: beginner 159776b2cf59SMatthew Knepley 159876b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 159976b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate() 160076b2cf59SMatthew Knepley @*/ 160176b2cf59SMatthew Knepley int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx) 160276b2cf59SMatthew Knepley { 160376b2cf59SMatthew Knepley PetscFunctionBegin; 160476b2cf59SMatthew Knepley PetscFunctionReturn(0); 160576b2cf59SMatthew Knepley } 160676b2cf59SMatthew Knepley 1607e74ef692SMatthew Knepley #undef __FUNCT__ 1608e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 1609ac226902SBarry Smith /*@C 161076b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 161176b2cf59SMatthew Knepley at the beginning of every step of the iteration. 161276b2cf59SMatthew Knepley 161376b2cf59SMatthew Knepley Collective on SNES 161476b2cf59SMatthew Knepley 161576b2cf59SMatthew Knepley Input Parameters: 161676b2cf59SMatthew Knepley . snes - The nonlinear solver context 161776b2cf59SMatthew Knepley . func - The function 161876b2cf59SMatthew Knepley 161976b2cf59SMatthew Knepley Calling sequence of func: 162076b2cf59SMatthew Knepley . func (TS ts, int step); 162176b2cf59SMatthew Knepley 162276b2cf59SMatthew Knepley . step - The current step of the iteration 162376b2cf59SMatthew Knepley 162476b2cf59SMatthew Knepley Level: intermediate 162576b2cf59SMatthew Knepley 162676b2cf59SMatthew Knepley .keywords: SNES, update 162776b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 162876b2cf59SMatthew Knepley @*/ 162976b2cf59SMatthew Knepley int SNESSetUpdate(SNES snes, int (*func)(SNES, int)) 163076b2cf59SMatthew Knepley { 163176b2cf59SMatthew Knepley PetscFunctionBegin; 16324482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 163376b2cf59SMatthew Knepley snes->update = func; 163476b2cf59SMatthew Knepley PetscFunctionReturn(0); 163576b2cf59SMatthew Knepley } 163676b2cf59SMatthew Knepley 1637e74ef692SMatthew Knepley #undef __FUNCT__ 1638e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 163976b2cf59SMatthew Knepley /*@ 164076b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 164176b2cf59SMatthew Knepley 164276b2cf59SMatthew Knepley Not collective 164376b2cf59SMatthew Knepley 164476b2cf59SMatthew Knepley Input Parameters: 164576b2cf59SMatthew Knepley . snes - The nonlinear solver context 164676b2cf59SMatthew Knepley . step - The current step of the iteration 164776b2cf59SMatthew Knepley 1648205452f4SMatthew Knepley Level: intermediate 1649205452f4SMatthew Knepley 165076b2cf59SMatthew Knepley .keywords: SNES, update 165176b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 165276b2cf59SMatthew Knepley @*/ 165376b2cf59SMatthew Knepley int SNESDefaultUpdate(SNES snes, int step) 165476b2cf59SMatthew Knepley { 165576b2cf59SMatthew Knepley PetscFunctionBegin; 165676b2cf59SMatthew Knepley PetscFunctionReturn(0); 165776b2cf59SMatthew Knepley } 165876b2cf59SMatthew Knepley 16594a2ae208SSatish Balay #undef __FUNCT__ 16604a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 16619b94acceSBarry Smith /* 16629b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 16639b94acceSBarry Smith positive parameter delta. 16649b94acceSBarry Smith 16659b94acceSBarry Smith Input Parameters: 1666c7afd0dbSLois Curfman McInnes + snes - the SNES context 16679b94acceSBarry Smith . y - approximate solution of linear system 16689b94acceSBarry Smith . fnorm - 2-norm of current function 1669c7afd0dbSLois Curfman McInnes - delta - trust region size 16709b94acceSBarry Smith 16719b94acceSBarry Smith Output Parameters: 1672c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 16739b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 16749b94acceSBarry Smith region, and exceeds zero otherwise. 1675c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 16769b94acceSBarry Smith 16779b94acceSBarry Smith Note: 16784b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 16799b94acceSBarry Smith is set to be the maximum allowable step size. 16809b94acceSBarry Smith 16819b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16829b94acceSBarry Smith */ 1683064f8208SBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 16849b94acceSBarry Smith { 1685064f8208SBarry Smith PetscReal nrm; 1686ea709b57SSatish Balay PetscScalar cnorm; 16873a40ed3dSBarry Smith int ierr; 16883a40ed3dSBarry Smith 16893a40ed3dSBarry Smith PetscFunctionBegin; 16904482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16914482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1692c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,2); 1693184914b5SBarry Smith 1694064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1695064f8208SBarry Smith if (nrm > *delta) { 1696064f8208SBarry Smith nrm = *delta/nrm; 1697064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1698064f8208SBarry Smith cnorm = nrm; 1699329f5518SBarry Smith ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 17009b94acceSBarry Smith *ynorm = *delta; 17019b94acceSBarry Smith } else { 17029b94acceSBarry Smith *gpnorm = 0.0; 1703064f8208SBarry Smith *ynorm = nrm; 17049b94acceSBarry Smith } 17053a40ed3dSBarry Smith PetscFunctionReturn(0); 17069b94acceSBarry Smith } 17079b94acceSBarry Smith 17085968eb51SBarry Smith static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero", 17095968eb51SBarry Smith "not currently used", 17105968eb51SBarry Smith "line search failed", 17115968eb51SBarry Smith "reach maximum number of iterations", 17125968eb51SBarry Smith "function norm became NaN (not a number)", 17135968eb51SBarry Smith "not currently used", 17145968eb51SBarry Smith "number of function computations exceeded", 17155968eb51SBarry Smith "not currently used", 17165968eb51SBarry Smith "still iterating", 17175968eb51SBarry Smith "not currently used", 17185968eb51SBarry Smith "absolute size of function norm", 17195968eb51SBarry Smith "relative decrease in function norm", 17205968eb51SBarry Smith "step size is small", 17215968eb51SBarry Smith "not currently used", 17225968eb51SBarry Smith "not currently used", 17235968eb51SBarry Smith "small trust region"}; 17245968eb51SBarry Smith 17254a2ae208SSatish Balay #undef __FUNCT__ 17264a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 17279b94acceSBarry Smith /*@ 17289b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 172963a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 17309b94acceSBarry Smith 1731c7afd0dbSLois Curfman McInnes Collective on SNES 1732c7afd0dbSLois Curfman McInnes 1733b2002411SLois Curfman McInnes Input Parameters: 1734c7afd0dbSLois Curfman McInnes + snes - the SNES context 1735c7afd0dbSLois Curfman McInnes - x - the solution vector 17369b94acceSBarry Smith 1737b2002411SLois Curfman McInnes Notes: 17388ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 17398ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 17408ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 17418ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 17428ddd3da0SLois Curfman McInnes 174336851e7fSLois Curfman McInnes Level: beginner 174436851e7fSLois Curfman McInnes 17459b94acceSBarry Smith .keywords: SNES, nonlinear, solve 17469b94acceSBarry Smith 1747*3ab0aad5SBarry Smith .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs() 17489b94acceSBarry Smith @*/ 1749c9780b6fSBarry Smith int SNESSolve(SNES snes,Vec x) 17509b94acceSBarry Smith { 1751f1af5d2fSBarry Smith int ierr; 1752f1af5d2fSBarry Smith PetscTruth flg; 1753052efed2SBarry Smith 17543a40ed3dSBarry Smith PetscFunctionBegin; 17554482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17564482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1757c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 175829bbc08cSBarry Smith if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 175974637425SBarry Smith 176082bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1761c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 1762758f92a0SBarry Smith if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 1763d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 176450ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1765c9780b6fSBarry Smith ierr = (*(snes)->solve)(snes);CHKERRQ(ierr); 1766d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1767b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 17688b179ff8SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1769da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1770da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 17715968eb51SBarry Smith if (snes->printreason) { 17725968eb51SBarry Smith if (snes->reason > 0) { 17735968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 17745968eb51SBarry Smith } else { 17755968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 17765968eb51SBarry Smith } 17775968eb51SBarry Smith } 17785968eb51SBarry Smith 17793a40ed3dSBarry Smith PetscFunctionReturn(0); 17809b94acceSBarry Smith } 17819b94acceSBarry Smith 17829b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 17839b94acceSBarry Smith 17844a2ae208SSatish Balay #undef __FUNCT__ 17854a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 178682bf6240SBarry Smith /*@C 17874b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 17889b94acceSBarry Smith 1789fee21e36SBarry Smith Collective on SNES 1790fee21e36SBarry Smith 1791c7afd0dbSLois Curfman McInnes Input Parameters: 1792c7afd0dbSLois Curfman McInnes + snes - the SNES context 1793454a90a3SBarry Smith - type - a known method 1794c7afd0dbSLois Curfman McInnes 1795c7afd0dbSLois Curfman McInnes Options Database Key: 1796454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1797c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1798ae12b187SLois Curfman McInnes 17999b94acceSBarry Smith Notes: 1800e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 18014b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1802c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 18034b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1804c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 18059b94acceSBarry Smith 1806ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1807ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1808ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1809ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1810ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1811ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1812ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1813ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1814ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1815b0a32e0cSBarry Smith appropriate method. 181636851e7fSLois Curfman McInnes 181736851e7fSLois Curfman McInnes Level: intermediate 1818a703fe33SLois Curfman McInnes 1819454a90a3SBarry Smith .keywords: SNES, set, type 1820435da068SBarry Smith 1821435da068SBarry Smith .seealso: SNESType, SNESCreate() 1822435da068SBarry Smith 18239b94acceSBarry Smith @*/ 18240e33f6ddSBarry Smith int SNESSetType(SNES snes,const SNESType type) 18259b94acceSBarry Smith { 18260f5bd95cSBarry Smith int ierr,(*r)(SNES); 18276831982aSBarry Smith PetscTruth match; 18283a40ed3dSBarry Smith 18293a40ed3dSBarry Smith PetscFunctionBegin; 18304482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18314482741eSBarry Smith PetscValidCharPointer(type,2); 183282bf6240SBarry Smith 18336831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 18340f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 183582bf6240SBarry Smith 183682bf6240SBarry Smith if (snes->setupcalled) { 1837e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 183882bf6240SBarry Smith snes->data = 0; 183982bf6240SBarry Smith } 18407d1a2b2bSBarry Smith 18419b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 184282bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 184382bf6240SBarry Smith 1844b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 184582bf6240SBarry Smith 184629bbc08cSBarry Smith if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type); 18471548b14aSBarry Smith 1848606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 184982bf6240SBarry Smith snes->data = 0; 18503a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1851454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 185282bf6240SBarry Smith 18533a40ed3dSBarry Smith PetscFunctionReturn(0); 18549b94acceSBarry Smith } 18559b94acceSBarry Smith 1856a847f771SSatish Balay 18579b94acceSBarry Smith /* --------------------------------------------------------------------- */ 18584a2ae208SSatish Balay #undef __FUNCT__ 18594a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 18609b94acceSBarry Smith /*@C 18619b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1862f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 18639b94acceSBarry Smith 1864fee21e36SBarry Smith Not Collective 1865fee21e36SBarry Smith 186636851e7fSLois Curfman McInnes Level: advanced 186736851e7fSLois Curfman McInnes 18689b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 18699b94acceSBarry Smith 18709b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 18719b94acceSBarry Smith @*/ 1872cf256101SBarry Smith int SNESRegisterDestroy(void) 18739b94acceSBarry Smith { 187482bf6240SBarry Smith int ierr; 187582bf6240SBarry Smith 18763a40ed3dSBarry Smith PetscFunctionBegin; 187782bf6240SBarry Smith if (SNESList) { 1878b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 187982bf6240SBarry Smith SNESList = 0; 18809b94acceSBarry Smith } 18814c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 18823a40ed3dSBarry Smith PetscFunctionReturn(0); 18839b94acceSBarry Smith } 18849b94acceSBarry Smith 18854a2ae208SSatish Balay #undef __FUNCT__ 18864a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 18879b94acceSBarry Smith /*@C 18889a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 18899b94acceSBarry Smith 1890c7afd0dbSLois Curfman McInnes Not Collective 1891c7afd0dbSLois Curfman McInnes 18929b94acceSBarry Smith Input Parameter: 18934b0e389bSBarry Smith . snes - nonlinear solver context 18949b94acceSBarry Smith 18959b94acceSBarry Smith Output Parameter: 18963a7fca6bSBarry Smith . type - SNES method (a character string) 18979b94acceSBarry Smith 189836851e7fSLois Curfman McInnes Level: intermediate 189936851e7fSLois Curfman McInnes 1900454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 19019b94acceSBarry Smith @*/ 1902454a90a3SBarry Smith int SNESGetType(SNES snes,SNESType *type) 19039b94acceSBarry Smith { 19043a40ed3dSBarry Smith PetscFunctionBegin; 19054482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19064482741eSBarry Smith PetscValidPointer(type,2); 1907454a90a3SBarry Smith *type = snes->type_name; 19083a40ed3dSBarry Smith PetscFunctionReturn(0); 19099b94acceSBarry Smith } 19109b94acceSBarry Smith 19114a2ae208SSatish Balay #undef __FUNCT__ 19124a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 19139b94acceSBarry Smith /*@C 19149b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 19159b94acceSBarry Smith stored. 19169b94acceSBarry Smith 1917c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1918c7afd0dbSLois Curfman McInnes 19199b94acceSBarry Smith Input Parameter: 19209b94acceSBarry Smith . snes - the SNES context 19219b94acceSBarry Smith 19229b94acceSBarry Smith Output Parameter: 19239b94acceSBarry Smith . x - the solution 19249b94acceSBarry Smith 192536851e7fSLois Curfman McInnes Level: advanced 192636851e7fSLois Curfman McInnes 19279b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 19289b94acceSBarry Smith 19294b27c08aSLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetSolutionUpdate() 19309b94acceSBarry Smith @*/ 19319b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 19329b94acceSBarry Smith { 19333a40ed3dSBarry Smith PetscFunctionBegin; 19344482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19354482741eSBarry Smith PetscValidPointer(x,2); 19369b94acceSBarry Smith *x = snes->vec_sol_always; 19373a40ed3dSBarry Smith PetscFunctionReturn(0); 19389b94acceSBarry Smith } 19399b94acceSBarry Smith 19404a2ae208SSatish Balay #undef __FUNCT__ 19414a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 19429b94acceSBarry Smith /*@C 19439b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 19449b94acceSBarry Smith stored. 19459b94acceSBarry Smith 1946c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1947c7afd0dbSLois Curfman McInnes 19489b94acceSBarry Smith Input Parameter: 19499b94acceSBarry Smith . snes - the SNES context 19509b94acceSBarry Smith 19519b94acceSBarry Smith Output Parameter: 19529b94acceSBarry Smith . x - the solution update 19539b94acceSBarry Smith 195436851e7fSLois Curfman McInnes Level: advanced 195536851e7fSLois Curfman McInnes 19569b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 19579b94acceSBarry Smith 19589b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 19599b94acceSBarry Smith @*/ 19609b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 19619b94acceSBarry Smith { 19623a40ed3dSBarry Smith PetscFunctionBegin; 19634482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19644482741eSBarry Smith PetscValidPointer(x,2); 19659b94acceSBarry Smith *x = snes->vec_sol_update_always; 19663a40ed3dSBarry Smith PetscFunctionReturn(0); 19679b94acceSBarry Smith } 19689b94acceSBarry Smith 19694a2ae208SSatish Balay #undef __FUNCT__ 19704a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 19719b94acceSBarry Smith /*@C 19723638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 19739b94acceSBarry Smith 1974c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1975c7afd0dbSLois Curfman McInnes 19769b94acceSBarry Smith Input Parameter: 19779b94acceSBarry Smith . snes - the SNES context 19789b94acceSBarry Smith 19799b94acceSBarry Smith Output Parameter: 19807bf4e008SBarry Smith + r - the function (or PETSC_NULL) 198100036973SBarry Smith . ctx - the function context (or PETSC_NULL) 198200036973SBarry Smith - func - the function (or PETSC_NULL) 19839b94acceSBarry Smith 198436851e7fSLois Curfman McInnes Level: advanced 198536851e7fSLois Curfman McInnes 1986a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 19879b94acceSBarry Smith 19884b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 19899b94acceSBarry Smith @*/ 199000036973SBarry Smith int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*)) 19919b94acceSBarry Smith { 19923a40ed3dSBarry Smith PetscFunctionBegin; 19934482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19947bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 19957bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 199600036973SBarry Smith if (func) *func = snes->computefunction; 19973a40ed3dSBarry Smith PetscFunctionReturn(0); 19989b94acceSBarry Smith } 19999b94acceSBarry Smith 20004a2ae208SSatish Balay #undef __FUNCT__ 20014a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 20023c7409f5SSatish Balay /*@C 20033c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 2004d850072dSLois Curfman McInnes SNES options in the database. 20053c7409f5SSatish Balay 2006fee21e36SBarry Smith Collective on SNES 2007fee21e36SBarry Smith 2008c7afd0dbSLois Curfman McInnes Input Parameter: 2009c7afd0dbSLois Curfman McInnes + snes - the SNES context 2010c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2011c7afd0dbSLois Curfman McInnes 2012d850072dSLois Curfman McInnes Notes: 2013a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2014c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2015d850072dSLois Curfman McInnes 201636851e7fSLois Curfman McInnes Level: advanced 201736851e7fSLois Curfman McInnes 20183c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 2019a86d99e1SLois Curfman McInnes 2020a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 20213c7409f5SSatish Balay @*/ 20221836bdbcSSatish Balay int SNESSetOptionsPrefix(SNES snes,const char prefix[]) 20233c7409f5SSatish Balay { 20243c7409f5SSatish Balay int ierr; 20253c7409f5SSatish Balay 20263a40ed3dSBarry Smith PetscFunctionBegin; 20274482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2028639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 202994b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20303a40ed3dSBarry Smith PetscFunctionReturn(0); 20313c7409f5SSatish Balay } 20323c7409f5SSatish Balay 20334a2ae208SSatish Balay #undef __FUNCT__ 20344a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 20353c7409f5SSatish Balay /*@C 2036f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2037d850072dSLois Curfman McInnes SNES options in the database. 20383c7409f5SSatish Balay 2039fee21e36SBarry Smith Collective on SNES 2040fee21e36SBarry Smith 2041c7afd0dbSLois Curfman McInnes Input Parameters: 2042c7afd0dbSLois Curfman McInnes + snes - the SNES context 2043c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2044c7afd0dbSLois Curfman McInnes 2045d850072dSLois Curfman McInnes Notes: 2046a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2047c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2048d850072dSLois Curfman McInnes 204936851e7fSLois Curfman McInnes Level: advanced 205036851e7fSLois Curfman McInnes 20513c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2052a86d99e1SLois Curfman McInnes 2053a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20543c7409f5SSatish Balay @*/ 20551836bdbcSSatish Balay int SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 20563c7409f5SSatish Balay { 20573c7409f5SSatish Balay int ierr; 20583c7409f5SSatish Balay 20593a40ed3dSBarry Smith PetscFunctionBegin; 20604482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2061639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 206294b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20633a40ed3dSBarry Smith PetscFunctionReturn(0); 20643c7409f5SSatish Balay } 20653c7409f5SSatish Balay 20664a2ae208SSatish Balay #undef __FUNCT__ 20674a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 20689ab63eb5SSatish Balay /*@C 20693c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20703c7409f5SSatish Balay SNES options in the database. 20713c7409f5SSatish Balay 2072c7afd0dbSLois Curfman McInnes Not Collective 2073c7afd0dbSLois Curfman McInnes 20743c7409f5SSatish Balay Input Parameter: 20753c7409f5SSatish Balay . snes - the SNES context 20763c7409f5SSatish Balay 20773c7409f5SSatish Balay Output Parameter: 20783c7409f5SSatish Balay . prefix - pointer to the prefix string used 20793c7409f5SSatish Balay 20809ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 20819ab63eb5SSatish Balay sufficient length to hold the prefix. 20829ab63eb5SSatish Balay 208336851e7fSLois Curfman McInnes Level: advanced 208436851e7fSLois Curfman McInnes 20853c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2086a86d99e1SLois Curfman McInnes 2087a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20883c7409f5SSatish Balay @*/ 20891836bdbcSSatish Balay int SNESGetOptionsPrefix(SNES snes,char *prefix[]) 20903c7409f5SSatish Balay { 20913c7409f5SSatish Balay int ierr; 20923c7409f5SSatish Balay 20933a40ed3dSBarry Smith PetscFunctionBegin; 20944482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2095639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 20963a40ed3dSBarry Smith PetscFunctionReturn(0); 20973c7409f5SSatish Balay } 20983c7409f5SSatish Balay 2099b2002411SLois Curfman McInnes 21004a2ae208SSatish Balay #undef __FUNCT__ 21014a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 21023cea93caSBarry Smith /*@C 21033cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 21043cea93caSBarry Smith 21057f6c08e0SMatthew Knepley Level: advanced 21063cea93caSBarry Smith @*/ 21071836bdbcSSatish Balay int SNESRegister(const char sname[],const char path[],const char name[],int (*function)(SNES)) 2108b2002411SLois Curfman McInnes { 2109b2002411SLois Curfman McInnes char fullname[256]; 2110b2002411SLois Curfman McInnes int ierr; 2111b2002411SLois Curfman McInnes 2112b2002411SLois Curfman McInnes PetscFunctionBegin; 2113b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2114c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2115b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2116b2002411SLois Curfman McInnes } 2117da9b6338SBarry Smith 2118da9b6338SBarry Smith #undef __FUNCT__ 2119da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 2120da9b6338SBarry Smith int SNESTestLocalMin(SNES snes) 2121da9b6338SBarry Smith { 2122da9b6338SBarry Smith int ierr,N,i,j; 2123da9b6338SBarry Smith Vec u,uh,fh; 2124da9b6338SBarry Smith PetscScalar value; 2125da9b6338SBarry Smith PetscReal norm; 2126da9b6338SBarry Smith 2127da9b6338SBarry Smith PetscFunctionBegin; 2128da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2129da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2130da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2131da9b6338SBarry Smith 2132da9b6338SBarry Smith /* currently only works for sequential */ 2133da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2134da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2135da9b6338SBarry Smith for (i=0; i<N; i++) { 2136da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 2137da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %d\n",i);CHKERRQ(ierr); 2138da9b6338SBarry Smith for (j=-10; j<11; j++) { 2139ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2140da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2141*3ab0aad5SBarry Smith ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2142da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2143da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %d %18.16e\n",j,norm);CHKERRQ(ierr); 2144da9b6338SBarry Smith value = -value; 2145da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2146da9b6338SBarry Smith } 2147da9b6338SBarry Smith } 2148da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2149da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2150da9b6338SBarry Smith PetscFunctionReturn(0); 2151da9b6338SBarry Smith } 2152