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 162*5968eb51SBarry Smith . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 163*5968eb51SBarry 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); 220*5968eb51SBarry Smith ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr); 221*5968eb51SBarry Smith if (flg) { 222*5968eb51SBarry Smith snes->printreason = PETSC_TRUE; 223*5968eb51SBarry 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 6649b94acceSBarry Smith /* --------------------------------------------------------------- */ 6654a2ae208SSatish Balay #undef __FUNCT__ 6664a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction" 6679b94acceSBarry Smith /*@C 6689b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6699b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6709b94acceSBarry Smith equations. 6719b94acceSBarry Smith 672fee21e36SBarry Smith Collective on SNES 673fee21e36SBarry Smith 674c7afd0dbSLois Curfman McInnes Input Parameters: 675c7afd0dbSLois Curfman McInnes + snes - the SNES context 676c7afd0dbSLois Curfman McInnes . func - function evaluation routine 677c7afd0dbSLois Curfman McInnes . r - vector to store function value 678c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 679c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6809b94acceSBarry Smith 681c7afd0dbSLois Curfman McInnes Calling sequence of func: 6828d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 683c7afd0dbSLois Curfman McInnes 684313e4042SLois Curfman McInnes . f - function vector 685c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 6869b94acceSBarry Smith 6879b94acceSBarry Smith Notes: 6889b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6899b94acceSBarry Smith $ f'(x) x = -f(x), 690c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 6919b94acceSBarry Smith 69236851e7fSLois Curfman McInnes Level: beginner 69336851e7fSLois Curfman McInnes 6949b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6959b94acceSBarry Smith 696a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6979b94acceSBarry Smith @*/ 6985334005bSBarry Smith int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 6999b94acceSBarry Smith { 7003a40ed3dSBarry Smith PetscFunctionBegin; 7014482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7024482741eSBarry Smith PetscValidHeaderSpecific(r,VEC_COOKIE,2); 703c9780b6fSBarry Smith PetscCheckSameComm(snes,1,r,2); 704184914b5SBarry Smith 7059b94acceSBarry Smith snes->computefunction = func; 7069b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7079b94acceSBarry Smith snes->funP = ctx; 7083a40ed3dSBarry Smith PetscFunctionReturn(0); 7099b94acceSBarry Smith } 7109b94acceSBarry Smith 7114a2ae208SSatish Balay #undef __FUNCT__ 7124a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction" 7139b94acceSBarry Smith /*@ 71436851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7159b94acceSBarry Smith SNESSetFunction(). 7169b94acceSBarry Smith 717c7afd0dbSLois Curfman McInnes Collective on SNES 718c7afd0dbSLois Curfman McInnes 7199b94acceSBarry Smith Input Parameters: 720c7afd0dbSLois Curfman McInnes + snes - the SNES context 721c7afd0dbSLois Curfman McInnes - x - input vector 7229b94acceSBarry Smith 7239b94acceSBarry Smith Output Parameter: 7243638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7259b94acceSBarry Smith 7261bffabb2SLois Curfman McInnes Notes: 72736851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 72836851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 72936851e7fSLois Curfman McInnes themselves. 73036851e7fSLois Curfman McInnes 73136851e7fSLois Curfman McInnes Level: developer 73236851e7fSLois Curfman McInnes 7339b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 7349b94acceSBarry Smith 735a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 7369b94acceSBarry Smith @*/ 7379b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x,Vec y) 7389b94acceSBarry Smith { 7399b94acceSBarry Smith int ierr; 7409b94acceSBarry Smith 7413a40ed3dSBarry Smith PetscFunctionBegin; 7424482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7434482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 7444482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,3); 745c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 746c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,3); 747184914b5SBarry Smith 748d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 749d64ed03dSBarry Smith PetscStackPush("SNES user function"); 7509b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 751d64ed03dSBarry Smith PetscStackPop; 752ae3c334cSLois Curfman McInnes snes->nfuncs++; 753d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 7543a40ed3dSBarry Smith PetscFunctionReturn(0); 7559b94acceSBarry Smith } 7569b94acceSBarry Smith 7574a2ae208SSatish Balay #undef __FUNCT__ 7584a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian" 75962fef451SLois Curfman McInnes /*@ 76062fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 76162fef451SLois Curfman McInnes set with SNESSetJacobian(). 76262fef451SLois Curfman McInnes 763c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 764c7afd0dbSLois Curfman McInnes 76562fef451SLois Curfman McInnes Input Parameters: 766c7afd0dbSLois Curfman McInnes + snes - the SNES context 767c7afd0dbSLois Curfman McInnes - x - input vector 76862fef451SLois Curfman McInnes 76962fef451SLois Curfman McInnes Output Parameters: 770c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 77162fef451SLois Curfman McInnes . B - optional preconditioning matrix 772c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 773fee21e36SBarry Smith 77462fef451SLois Curfman McInnes Notes: 77562fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 77662fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 77762fef451SLois Curfman McInnes 77894b7f48cSBarry Smith See KSPSetOperators() for important information about setting the 779dc5a77f8SLois Curfman McInnes flag parameter. 78062fef451SLois Curfman McInnes 78136851e7fSLois Curfman McInnes Level: developer 78236851e7fSLois Curfman McInnes 78362fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 78462fef451SLois Curfman McInnes 78594b7f48cSBarry Smith .seealso: SNESSetJacobian(), KSPSetOperators() 78662fef451SLois Curfman McInnes @*/ 7879b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 7889b94acceSBarry Smith { 7899b94acceSBarry Smith int ierr; 7903a40ed3dSBarry Smith 7913a40ed3dSBarry Smith PetscFunctionBegin; 7924482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 7934482741eSBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE,2); 7944482741eSBarry Smith PetscValidPointer(flg,5); 795c9780b6fSBarry Smith PetscCheckSameComm(snes,1,X,2); 7963a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 797d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 798c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 799d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 8009b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 801d64ed03dSBarry Smith PetscStackPop; 802d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 8036d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 8044482741eSBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 8054482741eSBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 8063a40ed3dSBarry Smith PetscFunctionReturn(0); 8079b94acceSBarry Smith } 8089b94acceSBarry Smith 8094a2ae208SSatish Balay #undef __FUNCT__ 8104a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian" 8119b94acceSBarry Smith /*@C 8129b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 813044dda88SLois Curfman McInnes location to store the matrix. 8149b94acceSBarry Smith 815c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 816c7afd0dbSLois Curfman McInnes 8179b94acceSBarry Smith Input Parameters: 818c7afd0dbSLois Curfman McInnes + snes - the SNES context 8199b94acceSBarry Smith . A - Jacobian matrix 8209b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 8219b94acceSBarry Smith . func - Jacobian evaluation routine 822c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 8232cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 8249b94acceSBarry Smith 8259b94acceSBarry Smith Calling sequence of func: 8268d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 8279b94acceSBarry Smith 828c7afd0dbSLois Curfman McInnes + x - input vector 8299b94acceSBarry Smith . A - Jacobian matrix 8309b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 831ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 83294b7f48cSBarry Smith structure (same as flag in KSPSetOperators()) 833c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 8349b94acceSBarry Smith 8359b94acceSBarry Smith Notes: 83694b7f48cSBarry Smith See KSPSetOperators() for important information about setting the flag 8372cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 838ac21db08SLois Curfman McInnes 839ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 8409b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 8419b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 8429b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 8439b94acceSBarry Smith throughout the global iterations. 8449b94acceSBarry Smith 84536851e7fSLois Curfman McInnes Level: beginner 84636851e7fSLois Curfman McInnes 8479b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 8489b94acceSBarry Smith 84994b7f48cSBarry Smith .seealso: KSPSetOperators(), SNESSetFunction() 8509b94acceSBarry Smith @*/ 851454a90a3SBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 8529b94acceSBarry Smith { 8533a7fca6bSBarry Smith int ierr; 8543a7fca6bSBarry Smith 8553a40ed3dSBarry Smith PetscFunctionBegin; 8564482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 8574482741eSBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 8584482741eSBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 859c9780b6fSBarry Smith if (A) PetscCheckSameComm(snes,1,A,2); 860c9780b6fSBarry Smith if (B) PetscCheckSameComm(snes,1,B,2); 8613a7fca6bSBarry Smith if (func) snes->computejacobian = func; 8623a7fca6bSBarry Smith if (ctx) snes->jacP = ctx; 8633a7fca6bSBarry Smith if (A) { 8643a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 8659b94acceSBarry Smith snes->jacobian = A; 8663a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 8673a7fca6bSBarry Smith } 8683a7fca6bSBarry Smith if (B) { 8693a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 8709b94acceSBarry Smith snes->jacobian_pre = B; 8713a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 8723a7fca6bSBarry Smith } 8733a40ed3dSBarry Smith PetscFunctionReturn(0); 8749b94acceSBarry Smith } 87562fef451SLois Curfman McInnes 8764a2ae208SSatish Balay #undef __FUNCT__ 8774a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian" 878c2aafc4cSSatish Balay /*@C 879b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 880b4fd4287SBarry Smith provided context for evaluating the Jacobian. 881b4fd4287SBarry Smith 882c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 883c7afd0dbSLois Curfman McInnes 884b4fd4287SBarry Smith Input Parameter: 885b4fd4287SBarry Smith . snes - the nonlinear solver context 886b4fd4287SBarry Smith 887b4fd4287SBarry Smith Output Parameters: 888c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 889b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 89000036973SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 89100036973SBarry Smith - func - location to put Jacobian function (or PETSC_NULL) 892fee21e36SBarry Smith 89336851e7fSLois Curfman McInnes Level: advanced 89436851e7fSLois Curfman McInnes 895b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 896b4fd4287SBarry Smith @*/ 89700036973SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 898b4fd4287SBarry Smith { 8993a40ed3dSBarry Smith PetscFunctionBegin; 9004482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 901b4fd4287SBarry Smith if (A) *A = snes->jacobian; 902b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 903b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 90400036973SBarry Smith if (func) *func = snes->computejacobian; 9053a40ed3dSBarry Smith PetscFunctionReturn(0); 906b4fd4287SBarry Smith } 907b4fd4287SBarry Smith 9089b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 90945fc7adcSBarry Smith extern int SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 9109b94acceSBarry Smith 9114a2ae208SSatish Balay #undef __FUNCT__ 9124a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp" 9139b94acceSBarry Smith /*@ 9149b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 915272ac6f2SLois Curfman McInnes of a nonlinear solver. 9169b94acceSBarry Smith 917fee21e36SBarry Smith Collective on SNES 918fee21e36SBarry Smith 919c7afd0dbSLois Curfman McInnes Input Parameters: 920c7afd0dbSLois Curfman McInnes + snes - the SNES context 921c7afd0dbSLois Curfman McInnes - x - the solution vector 922c7afd0dbSLois Curfman McInnes 923272ac6f2SLois Curfman McInnes Notes: 924272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 925272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 926272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 927272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 928272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 929272ac6f2SLois Curfman McInnes 93036851e7fSLois Curfman McInnes Level: advanced 93136851e7fSLois Curfman McInnes 9329b94acceSBarry Smith .keywords: SNES, nonlinear, setup 9339b94acceSBarry Smith 9349b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 9359b94acceSBarry Smith @*/ 9368ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 9379b94acceSBarry Smith { 938f1af5d2fSBarry Smith int ierr; 9394b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 9403a40ed3dSBarry Smith 9413a40ed3dSBarry Smith PetscFunctionBegin; 9424482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 9434482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 944c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 9458ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 9469b94acceSBarry Smith 947b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 948c1f60f51SBarry Smith /* 949c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 950dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 951c1f60f51SBarry Smith */ 952c1f60f51SBarry Smith if (flg) { 953c1f60f51SBarry Smith Mat J; 9545a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 9555a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 9563a7fca6bSBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n"); 9573a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 9583a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 959c1f60f51SBarry Smith } 96045fc7adcSBarry Smith 961b4014bb2SMatthew Knepley #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) 96245fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 96345fc7adcSBarry Smith if (flg) { 96445fc7adcSBarry Smith Mat J; 96545fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 96645fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 96745fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 96845fc7adcSBarry Smith } 96932a4b47aSMatthew Knepley #endif 97045fc7adcSBarry Smith 971b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 972c1f60f51SBarry Smith /* 973dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 974c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 975c1f60f51SBarry Smith */ 97631615d04SLois Curfman McInnes if (flg) { 977272ac6f2SLois Curfman McInnes Mat J; 978b5d62d44SBarry Smith KSP ksp; 97994b7f48cSBarry Smith PC pc; 980f3ef73ceSBarry Smith 9815a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 9823a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 98393b98531SBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"); 9843a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 9853a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 9863a7fca6bSBarry Smith 987f3ef73ceSBarry Smith /* force no preconditioner */ 98894b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 989b5d62d44SBarry Smith ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 990a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 991a9815358SBarry Smith if (!flg) { 992f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 993272ac6f2SLois Curfman McInnes } 994a9815358SBarry Smith } 995f3ef73ceSBarry Smith 99629bbc08cSBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 99729bbc08cSBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 99829bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 999a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 100029bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1001a8c6a408SBarry Smith } 1002a703fe33SLois Curfman McInnes 1003a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 10044b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 10056831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 100694b7f48cSBarry Smith KSP ksp; 100794b7f48cSBarry Smith ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1008186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1009a703fe33SLois Curfman McInnes } 10104b27c08aSLois Curfman McInnes 1011a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 101282bf6240SBarry Smith snes->setupcalled = 1; 10133a40ed3dSBarry Smith PetscFunctionReturn(0); 10149b94acceSBarry Smith } 10159b94acceSBarry Smith 10164a2ae208SSatish Balay #undef __FUNCT__ 10174a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 10189b94acceSBarry Smith /*@C 10199b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 10209b94acceSBarry Smith with SNESCreate(). 10219b94acceSBarry Smith 1022c7afd0dbSLois Curfman McInnes Collective on SNES 1023c7afd0dbSLois Curfman McInnes 10249b94acceSBarry Smith Input Parameter: 10259b94acceSBarry Smith . snes - the SNES context 10269b94acceSBarry Smith 102736851e7fSLois Curfman McInnes Level: beginner 102836851e7fSLois Curfman McInnes 10299b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 10309b94acceSBarry Smith 103163a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 10329b94acceSBarry Smith @*/ 10339b94acceSBarry Smith int SNESDestroy(SNES snes) 10349b94acceSBarry Smith { 1035b8a78c4aSBarry Smith int i,ierr; 10363a40ed3dSBarry Smith 10373a40ed3dSBarry Smith PetscFunctionBegin; 10384482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 10393a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1040d4bb536fSBarry Smith 1041be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 10420f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1043be0abb6dSBarry Smith 1044e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1045606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 10463a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 10473a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 104894b7f48cSBarry Smith ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1049522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1050b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1051b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1052b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1053b8a78c4aSBarry Smith } 1054b8a78c4aSBarry Smith } 1055b0a32e0cSBarry Smith PetscLogObjectDestroy((PetscObject)snes); 10560452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 10573a40ed3dSBarry Smith PetscFunctionReturn(0); 10589b94acceSBarry Smith } 10599b94acceSBarry Smith 10609b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 10619b94acceSBarry Smith 10624a2ae208SSatish Balay #undef __FUNCT__ 10634a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 10649b94acceSBarry Smith /*@ 1065d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 10669b94acceSBarry Smith 1067c7afd0dbSLois Curfman McInnes Collective on SNES 1068c7afd0dbSLois Curfman McInnes 10699b94acceSBarry Smith Input Parameters: 1070c7afd0dbSLois Curfman McInnes + snes - the SNES context 107133174efeSLois Curfman McInnes . atol - absolute convergence tolerance 107233174efeSLois Curfman McInnes . rtol - relative convergence tolerance 107333174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 107433174efeSLois Curfman McInnes of the change in the solution between steps 107533174efeSLois Curfman McInnes . maxit - maximum number of iterations 1076c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1077fee21e36SBarry Smith 107833174efeSLois Curfman McInnes Options Database Keys: 1079c7afd0dbSLois Curfman McInnes + -snes_atol <atol> - Sets atol 1080c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1081c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1082c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1083c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 10849b94acceSBarry Smith 1085d7a720efSLois Curfman McInnes Notes: 10869b94acceSBarry Smith The default maximum number of iterations is 50. 10879b94acceSBarry Smith The default maximum number of function evaluations is 1000. 10889b94acceSBarry Smith 108936851e7fSLois Curfman McInnes Level: intermediate 109036851e7fSLois Curfman McInnes 109133174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 10929b94acceSBarry Smith 1093d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 10949b94acceSBarry Smith @*/ 1095329f5518SBarry Smith int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf) 10969b94acceSBarry Smith { 10973a40ed3dSBarry Smith PetscFunctionBegin; 10984482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1099d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1100d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1101d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1102d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1103d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11043a40ed3dSBarry Smith PetscFunctionReturn(0); 11059b94acceSBarry Smith } 11069b94acceSBarry Smith 11074a2ae208SSatish Balay #undef __FUNCT__ 11084a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11099b94acceSBarry Smith /*@ 111033174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 111133174efeSLois Curfman McInnes 1112c7afd0dbSLois Curfman McInnes Not Collective 1113c7afd0dbSLois Curfman McInnes 111433174efeSLois Curfman McInnes Input Parameters: 1115c7afd0dbSLois Curfman McInnes + snes - the SNES context 111633174efeSLois Curfman McInnes . atol - absolute convergence tolerance 111733174efeSLois Curfman McInnes . rtol - relative convergence tolerance 111833174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 111933174efeSLois Curfman McInnes of the change in the solution between steps 112033174efeSLois Curfman McInnes . maxit - maximum number of iterations 1121c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1122fee21e36SBarry Smith 112333174efeSLois Curfman McInnes Notes: 112433174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 112533174efeSLois Curfman McInnes 112636851e7fSLois Curfman McInnes Level: intermediate 112736851e7fSLois Curfman McInnes 112833174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 112933174efeSLois Curfman McInnes 113033174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 113133174efeSLois Curfman McInnes @*/ 1132329f5518SBarry Smith int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf) 113333174efeSLois Curfman McInnes { 11343a40ed3dSBarry Smith PetscFunctionBegin; 11354482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 113633174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 113733174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 113833174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 113933174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 114033174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 11413a40ed3dSBarry Smith PetscFunctionReturn(0); 114233174efeSLois Curfman McInnes } 114333174efeSLois Curfman McInnes 11444a2ae208SSatish Balay #undef __FUNCT__ 11454a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 114633174efeSLois Curfman McInnes /*@ 11479b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 11489b94acceSBarry Smith 1149fee21e36SBarry Smith Collective on SNES 1150fee21e36SBarry Smith 1151c7afd0dbSLois Curfman McInnes Input Parameters: 1152c7afd0dbSLois Curfman McInnes + snes - the SNES context 1153c7afd0dbSLois Curfman McInnes - tol - tolerance 1154c7afd0dbSLois Curfman McInnes 11559b94acceSBarry Smith Options Database Key: 1156c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 11579b94acceSBarry Smith 115836851e7fSLois Curfman McInnes Level: intermediate 115936851e7fSLois Curfman McInnes 11609b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 11619b94acceSBarry Smith 1162d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 11639b94acceSBarry Smith @*/ 1164329f5518SBarry Smith int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 11659b94acceSBarry Smith { 11663a40ed3dSBarry Smith PetscFunctionBegin; 11674482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 11689b94acceSBarry Smith snes->deltatol = tol; 11693a40ed3dSBarry Smith PetscFunctionReturn(0); 11709b94acceSBarry Smith } 11719b94acceSBarry Smith 1172df9fa365SBarry Smith /* 1173df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1174df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1175df9fa365SBarry Smith macros instead of functions 1176df9fa365SBarry Smith */ 11774a2ae208SSatish Balay #undef __FUNCT__ 11784a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 1179329f5518SBarry Smith int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx) 1180ce1608b8SBarry Smith { 1181ce1608b8SBarry Smith int ierr; 1182ce1608b8SBarry Smith 1183ce1608b8SBarry Smith PetscFunctionBegin; 11844482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1185ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1186ce1608b8SBarry Smith PetscFunctionReturn(0); 1187ce1608b8SBarry Smith } 1188ce1608b8SBarry Smith 11894a2ae208SSatish Balay #undef __FUNCT__ 11904a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 11911836bdbcSSatish Balay int SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1192df9fa365SBarry Smith { 1193df9fa365SBarry Smith int ierr; 1194df9fa365SBarry Smith 1195df9fa365SBarry Smith PetscFunctionBegin; 1196df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1197df9fa365SBarry Smith PetscFunctionReturn(0); 1198df9fa365SBarry Smith } 1199df9fa365SBarry Smith 12004a2ae208SSatish Balay #undef __FUNCT__ 12014a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 1202b0a32e0cSBarry Smith int SNESLGMonitorDestroy(PetscDrawLG draw) 1203df9fa365SBarry Smith { 1204df9fa365SBarry Smith int ierr; 1205df9fa365SBarry Smith 1206df9fa365SBarry Smith PetscFunctionBegin; 1207df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1208df9fa365SBarry Smith PetscFunctionReturn(0); 1209df9fa365SBarry Smith } 1210df9fa365SBarry Smith 12119b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12129b94acceSBarry Smith 12134a2ae208SSatish Balay #undef __FUNCT__ 12144a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12159b94acceSBarry Smith /*@C 12165cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12179b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12189b94acceSBarry Smith progress. 12199b94acceSBarry Smith 1220fee21e36SBarry Smith Collective on SNES 1221fee21e36SBarry Smith 1222c7afd0dbSLois Curfman McInnes Input Parameters: 1223c7afd0dbSLois Curfman McInnes + snes - the SNES context 1224c7afd0dbSLois Curfman McInnes . func - monitoring routine 1225b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1226b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desitre) 1227b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1228b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 12299b94acceSBarry Smith 1230c7afd0dbSLois Curfman McInnes Calling sequence of func: 1231329f5518SBarry Smith $ int func(SNES snes,int its, PetscReal norm,void *mctx) 1232c7afd0dbSLois Curfman McInnes 1233c7afd0dbSLois Curfman McInnes + snes - the SNES context 1234c7afd0dbSLois Curfman McInnes . its - iteration number 1235c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 123640a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 12379b94acceSBarry Smith 12389665c990SLois Curfman McInnes Options Database Keys: 1239c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1240c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1241c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1242c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1243c7afd0dbSLois Curfman McInnes been hardwired into a code by 1244c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1245c7afd0dbSLois Curfman McInnes does not cancel those set via 1246c7afd0dbSLois Curfman McInnes the options database. 12479665c990SLois Curfman McInnes 1248639f9d9dSBarry Smith Notes: 12496bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 12506bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 12516bc08f3fSLois Curfman McInnes order in which they were set. 1252639f9d9dSBarry Smith 125336851e7fSLois Curfman McInnes Level: intermediate 125436851e7fSLois Curfman McInnes 12559b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 12569b94acceSBarry Smith 12575cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 12589b94acceSBarry Smith @*/ 1259329f5518SBarry Smith int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *)) 12609b94acceSBarry Smith { 12613a40ed3dSBarry Smith PetscFunctionBegin; 12624482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1263639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 126429bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1265639f9d9dSBarry Smith } 1266639f9d9dSBarry Smith 1267639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1268b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1269639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 12703a40ed3dSBarry Smith PetscFunctionReturn(0); 12719b94acceSBarry Smith } 12729b94acceSBarry Smith 12734a2ae208SSatish Balay #undef __FUNCT__ 12744a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 12755cd90555SBarry Smith /*@C 12765cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 12775cd90555SBarry Smith 1278c7afd0dbSLois Curfman McInnes Collective on SNES 1279c7afd0dbSLois Curfman McInnes 12805cd90555SBarry Smith Input Parameters: 12815cd90555SBarry Smith . snes - the SNES context 12825cd90555SBarry Smith 12831a480d89SAdministrator Options Database Key: 1284c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1285c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1286c7afd0dbSLois Curfman McInnes set via the options database 12875cd90555SBarry Smith 12885cd90555SBarry Smith Notes: 12895cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 12905cd90555SBarry Smith 129136851e7fSLois Curfman McInnes Level: intermediate 129236851e7fSLois Curfman McInnes 12935cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 12945cd90555SBarry Smith 12955cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 12965cd90555SBarry Smith @*/ 12975cd90555SBarry Smith int SNESClearMonitor(SNES snes) 12985cd90555SBarry Smith { 12995cd90555SBarry Smith PetscFunctionBegin; 13004482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13015cd90555SBarry Smith snes->numbermonitors = 0; 13025cd90555SBarry Smith PetscFunctionReturn(0); 13035cd90555SBarry Smith } 13045cd90555SBarry Smith 13054a2ae208SSatish Balay #undef __FUNCT__ 13064a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13079b94acceSBarry Smith /*@C 13089b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13099b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13109b94acceSBarry Smith 1311fee21e36SBarry Smith Collective on SNES 1312fee21e36SBarry Smith 1313c7afd0dbSLois Curfman McInnes Input Parameters: 1314c7afd0dbSLois Curfman McInnes + snes - the SNES context 1315c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1316c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1317c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 13189b94acceSBarry Smith 1319c7afd0dbSLois Curfman McInnes Calling sequence of func: 1320329f5518SBarry Smith $ int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1321c7afd0dbSLois Curfman McInnes 1322c7afd0dbSLois Curfman McInnes + snes - the SNES context 1323c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1324184914b5SBarry Smith . reason - reason for convergence/divergence 1325c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 13264b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 13274b27c08aSLois Curfman McInnes - f - 2-norm of function 13289b94acceSBarry Smith 132936851e7fSLois Curfman McInnes Level: advanced 133036851e7fSLois Curfman McInnes 13319b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13329b94acceSBarry Smith 13334b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 13349b94acceSBarry Smith @*/ 1335329f5518SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 13369b94acceSBarry Smith { 13373a40ed3dSBarry Smith PetscFunctionBegin; 13384482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13399b94acceSBarry Smith (snes)->converged = func; 13409b94acceSBarry Smith (snes)->cnvP = cctx; 13413a40ed3dSBarry Smith PetscFunctionReturn(0); 13429b94acceSBarry Smith } 13439b94acceSBarry Smith 13444a2ae208SSatish Balay #undef __FUNCT__ 13454a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 1346184914b5SBarry Smith /*@C 1347184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1348184914b5SBarry Smith 1349184914b5SBarry Smith Not Collective 1350184914b5SBarry Smith 1351184914b5SBarry Smith Input Parameter: 1352184914b5SBarry Smith . snes - the SNES context 1353184914b5SBarry Smith 1354184914b5SBarry Smith Output Parameter: 1355e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1356184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1357184914b5SBarry Smith 1358184914b5SBarry Smith Level: intermediate 1359184914b5SBarry Smith 1360184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1361184914b5SBarry Smith 1362184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1363184914b5SBarry Smith 13644b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1365184914b5SBarry Smith @*/ 1366184914b5SBarry Smith int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1367184914b5SBarry Smith { 1368184914b5SBarry Smith PetscFunctionBegin; 13694482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 13704482741eSBarry Smith PetscValidPointer(reason,2); 1371184914b5SBarry Smith *reason = snes->reason; 1372184914b5SBarry Smith PetscFunctionReturn(0); 1373184914b5SBarry Smith } 1374184914b5SBarry Smith 13754a2ae208SSatish Balay #undef __FUNCT__ 13764a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1377c9005455SLois Curfman McInnes /*@ 1378c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1379c9005455SLois Curfman McInnes 1380fee21e36SBarry Smith Collective on SNES 1381fee21e36SBarry Smith 1382c7afd0dbSLois Curfman McInnes Input Parameters: 1383c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1384c7afd0dbSLois Curfman McInnes . a - array to hold history 1385cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1386758f92a0SBarry Smith . na - size of a and its 138764731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1388758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1389c7afd0dbSLois Curfman McInnes 1390c9005455SLois Curfman McInnes Notes: 13914b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1392c9005455SLois Curfman McInnes at each step. 1393c9005455SLois Curfman McInnes 1394c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1395c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1396c9005455SLois Curfman McInnes during the section of code that is being timed. 1397c9005455SLois Curfman McInnes 139836851e7fSLois Curfman McInnes Level: intermediate 139936851e7fSLois Curfman McInnes 1400c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1401758f92a0SBarry Smith 140208405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1403758f92a0SBarry Smith 1404c9005455SLois Curfman McInnes @*/ 14051836bdbcSSatish Balay int SNESSetConvergenceHistory(SNES snes,PetscReal a[],int *its,int na,PetscTruth reset) 1406c9005455SLois Curfman McInnes { 14073a40ed3dSBarry Smith PetscFunctionBegin; 14084482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 14094482741eSBarry Smith if (na) PetscValidScalarPointer(a,2); 1410c9005455SLois Curfman McInnes snes->conv_hist = a; 1411758f92a0SBarry Smith snes->conv_hist_its = its; 1412758f92a0SBarry Smith snes->conv_hist_max = na; 1413758f92a0SBarry Smith snes->conv_hist_reset = reset; 1414758f92a0SBarry Smith PetscFunctionReturn(0); 1415758f92a0SBarry Smith } 1416758f92a0SBarry Smith 14174a2ae208SSatish Balay #undef __FUNCT__ 14184a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 14190c4c9dddSBarry Smith /*@C 1420758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1421758f92a0SBarry Smith 1422758f92a0SBarry Smith Collective on SNES 1423758f92a0SBarry Smith 1424758f92a0SBarry Smith Input Parameter: 1425758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1426758f92a0SBarry Smith 1427758f92a0SBarry Smith Output Parameters: 1428758f92a0SBarry Smith . a - array to hold history 1429758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1430758f92a0SBarry Smith negative if not converged) for each solve. 1431758f92a0SBarry Smith - na - size of a and its 1432758f92a0SBarry Smith 1433758f92a0SBarry Smith Notes: 1434758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1435758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1436758f92a0SBarry Smith 1437758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1438758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1439758f92a0SBarry Smith during the section of code that is being timed. 1440758f92a0SBarry Smith 1441758f92a0SBarry Smith Level: intermediate 1442758f92a0SBarry Smith 1443758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1444758f92a0SBarry Smith 1445758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1446758f92a0SBarry Smith 1447758f92a0SBarry Smith @*/ 14481836bdbcSSatish Balay int SNESGetConvergenceHistory(SNES snes,PetscReal *a[],int *its[],int *na) 1449758f92a0SBarry Smith { 1450758f92a0SBarry Smith PetscFunctionBegin; 14514482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1452758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1453758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1454758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 14553a40ed3dSBarry Smith PetscFunctionReturn(0); 1456c9005455SLois Curfman McInnes } 1457c9005455SLois Curfman McInnes 1458e74ef692SMatthew Knepley #undef __FUNCT__ 1459e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC" 1460ac226902SBarry Smith /*@C 146176b2cf59SMatthew Knepley SNESSetRhsBC - Sets the function which applies boundary conditions 146276b2cf59SMatthew Knepley to the Rhs of each system. 146376b2cf59SMatthew Knepley 146476b2cf59SMatthew Knepley Collective on SNES 146576b2cf59SMatthew Knepley 146676b2cf59SMatthew Knepley Input Parameters: 146776b2cf59SMatthew Knepley . snes - The nonlinear solver context 146876b2cf59SMatthew Knepley . func - The function 146976b2cf59SMatthew Knepley 147076b2cf59SMatthew Knepley Calling sequence of func: 147176b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx); 147276b2cf59SMatthew Knepley 147376b2cf59SMatthew Knepley . rhs - The current rhs vector 147476b2cf59SMatthew Knepley . ctx - The user-context 147576b2cf59SMatthew Knepley 147676b2cf59SMatthew Knepley Level: intermediate 147776b2cf59SMatthew Knepley 147876b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 147976b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate() 148076b2cf59SMatthew Knepley @*/ 148176b2cf59SMatthew Knepley int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *)) 148276b2cf59SMatthew Knepley { 148376b2cf59SMatthew Knepley PetscFunctionBegin; 14844482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 148576b2cf59SMatthew Knepley snes->applyrhsbc = func; 148676b2cf59SMatthew Knepley PetscFunctionReturn(0); 148776b2cf59SMatthew Knepley } 148876b2cf59SMatthew Knepley 1489e74ef692SMatthew Knepley #undef __FUNCT__ 1490e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC" 149176b2cf59SMatthew Knepley /*@ 149276b2cf59SMatthew Knepley SNESDefaultRhsBC - The default boundary condition function which does nothing. 149376b2cf59SMatthew Knepley 149476b2cf59SMatthew Knepley Not collective 149576b2cf59SMatthew Knepley 149676b2cf59SMatthew Knepley Input Parameters: 149776b2cf59SMatthew Knepley . snes - The nonlinear solver context 149876b2cf59SMatthew Knepley . rhs - The Rhs 149976b2cf59SMatthew Knepley . ctx - The user-context 150076b2cf59SMatthew Knepley 150176b2cf59SMatthew Knepley Level: beginner 150276b2cf59SMatthew Knepley 150376b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 150476b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate() 150576b2cf59SMatthew Knepley @*/ 150676b2cf59SMatthew Knepley int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx) 150776b2cf59SMatthew Knepley { 150876b2cf59SMatthew Knepley PetscFunctionBegin; 150976b2cf59SMatthew Knepley PetscFunctionReturn(0); 151076b2cf59SMatthew Knepley } 151176b2cf59SMatthew Knepley 1512e74ef692SMatthew Knepley #undef __FUNCT__ 1513e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC" 1514ac226902SBarry Smith /*@C 151576b2cf59SMatthew Knepley SNESSetSolutionBC - Sets the function which applies boundary conditions 151676b2cf59SMatthew Knepley to the solution of each system. 151776b2cf59SMatthew Knepley 151876b2cf59SMatthew Knepley Collective on SNES 151976b2cf59SMatthew Knepley 152076b2cf59SMatthew Knepley Input Parameters: 152176b2cf59SMatthew Knepley . snes - The nonlinear solver context 152276b2cf59SMatthew Knepley . func - The function 152376b2cf59SMatthew Knepley 152476b2cf59SMatthew Knepley Calling sequence of func: 15259d06fb6bSMatthew Knepley . func (SNES snes, Vec rsol, void *ctx); 152676b2cf59SMatthew Knepley 152776b2cf59SMatthew Knepley . sol - The current solution vector 152876b2cf59SMatthew Knepley . ctx - The user-context 152976b2cf59SMatthew Knepley 153076b2cf59SMatthew Knepley Level: intermediate 153176b2cf59SMatthew Knepley 153276b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 153376b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate() 153476b2cf59SMatthew Knepley @*/ 153576b2cf59SMatthew Knepley int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *)) 153676b2cf59SMatthew Knepley { 153776b2cf59SMatthew Knepley PetscFunctionBegin; 15384482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 153976b2cf59SMatthew Knepley snes->applysolbc = func; 154076b2cf59SMatthew Knepley PetscFunctionReturn(0); 154176b2cf59SMatthew Knepley } 154276b2cf59SMatthew Knepley 1543e74ef692SMatthew Knepley #undef __FUNCT__ 1544e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC" 154576b2cf59SMatthew Knepley /*@ 154676b2cf59SMatthew Knepley SNESDefaultSolutionBC - The default boundary condition function which does nothing. 154776b2cf59SMatthew Knepley 154876b2cf59SMatthew Knepley Not collective 154976b2cf59SMatthew Knepley 155076b2cf59SMatthew Knepley Input Parameters: 155176b2cf59SMatthew Knepley . snes - The nonlinear solver context 155276b2cf59SMatthew Knepley . sol - The solution 155376b2cf59SMatthew Knepley . ctx - The user-context 155476b2cf59SMatthew Knepley 155576b2cf59SMatthew Knepley Level: beginner 155676b2cf59SMatthew Knepley 155776b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 155876b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate() 155976b2cf59SMatthew Knepley @*/ 156076b2cf59SMatthew Knepley int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx) 156176b2cf59SMatthew Knepley { 156276b2cf59SMatthew Knepley PetscFunctionBegin; 156376b2cf59SMatthew Knepley PetscFunctionReturn(0); 156476b2cf59SMatthew Knepley } 156576b2cf59SMatthew Knepley 1566e74ef692SMatthew Knepley #undef __FUNCT__ 1567e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 1568ac226902SBarry Smith /*@C 156976b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 157076b2cf59SMatthew Knepley at the beginning of every step of the iteration. 157176b2cf59SMatthew Knepley 157276b2cf59SMatthew Knepley Collective on SNES 157376b2cf59SMatthew Knepley 157476b2cf59SMatthew Knepley Input Parameters: 157576b2cf59SMatthew Knepley . snes - The nonlinear solver context 157676b2cf59SMatthew Knepley . func - The function 157776b2cf59SMatthew Knepley 157876b2cf59SMatthew Knepley Calling sequence of func: 157976b2cf59SMatthew Knepley . func (TS ts, int step); 158076b2cf59SMatthew Knepley 158176b2cf59SMatthew Knepley . step - The current step of the iteration 158276b2cf59SMatthew Knepley 158376b2cf59SMatthew Knepley Level: intermediate 158476b2cf59SMatthew Knepley 158576b2cf59SMatthew Knepley .keywords: SNES, update 158676b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 158776b2cf59SMatthew Knepley @*/ 158876b2cf59SMatthew Knepley int SNESSetUpdate(SNES snes, int (*func)(SNES, int)) 158976b2cf59SMatthew Knepley { 159076b2cf59SMatthew Knepley PetscFunctionBegin; 15914482741eSBarry Smith PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 159276b2cf59SMatthew Knepley snes->update = func; 159376b2cf59SMatthew Knepley PetscFunctionReturn(0); 159476b2cf59SMatthew Knepley } 159576b2cf59SMatthew Knepley 1596e74ef692SMatthew Knepley #undef __FUNCT__ 1597e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 159876b2cf59SMatthew Knepley /*@ 159976b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 160076b2cf59SMatthew Knepley 160176b2cf59SMatthew Knepley Not collective 160276b2cf59SMatthew Knepley 160376b2cf59SMatthew Knepley Input Parameters: 160476b2cf59SMatthew Knepley . snes - The nonlinear solver context 160576b2cf59SMatthew Knepley . step - The current step of the iteration 160676b2cf59SMatthew Knepley 1607205452f4SMatthew Knepley Level: intermediate 1608205452f4SMatthew Knepley 160976b2cf59SMatthew Knepley .keywords: SNES, update 161076b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 161176b2cf59SMatthew Knepley @*/ 161276b2cf59SMatthew Knepley int SNESDefaultUpdate(SNES snes, int step) 161376b2cf59SMatthew Knepley { 161476b2cf59SMatthew Knepley PetscFunctionBegin; 161576b2cf59SMatthew Knepley PetscFunctionReturn(0); 161676b2cf59SMatthew Knepley } 161776b2cf59SMatthew Knepley 16184a2ae208SSatish Balay #undef __FUNCT__ 16194a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 16209b94acceSBarry Smith /* 16219b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 16229b94acceSBarry Smith positive parameter delta. 16239b94acceSBarry Smith 16249b94acceSBarry Smith Input Parameters: 1625c7afd0dbSLois Curfman McInnes + snes - the SNES context 16269b94acceSBarry Smith . y - approximate solution of linear system 16279b94acceSBarry Smith . fnorm - 2-norm of current function 1628c7afd0dbSLois Curfman McInnes - delta - trust region size 16299b94acceSBarry Smith 16309b94acceSBarry Smith Output Parameters: 1631c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 16329b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 16339b94acceSBarry Smith region, and exceeds zero otherwise. 1634c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 16359b94acceSBarry Smith 16369b94acceSBarry Smith Note: 16374b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 16389b94acceSBarry Smith is set to be the maximum allowable step size. 16399b94acceSBarry Smith 16409b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16419b94acceSBarry Smith */ 1642064f8208SBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 16439b94acceSBarry Smith { 1644064f8208SBarry Smith PetscReal nrm; 1645ea709b57SSatish Balay PetscScalar cnorm; 16463a40ed3dSBarry Smith int ierr; 16473a40ed3dSBarry Smith 16483a40ed3dSBarry Smith PetscFunctionBegin; 16494482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 16504482741eSBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1651c9780b6fSBarry Smith PetscCheckSameComm(snes,1,y,2); 1652184914b5SBarry Smith 1653064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1654064f8208SBarry Smith if (nrm > *delta) { 1655064f8208SBarry Smith nrm = *delta/nrm; 1656064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1657064f8208SBarry Smith cnorm = nrm; 1658329f5518SBarry Smith ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 16599b94acceSBarry Smith *ynorm = *delta; 16609b94acceSBarry Smith } else { 16619b94acceSBarry Smith *gpnorm = 0.0; 1662064f8208SBarry Smith *ynorm = nrm; 16639b94acceSBarry Smith } 16643a40ed3dSBarry Smith PetscFunctionReturn(0); 16659b94acceSBarry Smith } 16669b94acceSBarry Smith 1667*5968eb51SBarry Smith static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero", 1668*5968eb51SBarry Smith "not currently used", 1669*5968eb51SBarry Smith "line search failed", 1670*5968eb51SBarry Smith "reach maximum number of iterations", 1671*5968eb51SBarry Smith "function norm became NaN (not a number)", 1672*5968eb51SBarry Smith "not currently used", 1673*5968eb51SBarry Smith "number of function computations exceeded", 1674*5968eb51SBarry Smith "not currently used", 1675*5968eb51SBarry Smith "still iterating", 1676*5968eb51SBarry Smith "not currently used", 1677*5968eb51SBarry Smith "absolute size of function norm", 1678*5968eb51SBarry Smith "relative decrease in function norm", 1679*5968eb51SBarry Smith "step size is small", 1680*5968eb51SBarry Smith "not currently used", 1681*5968eb51SBarry Smith "not currently used", 1682*5968eb51SBarry Smith "small trust region"}; 1683*5968eb51SBarry Smith 16844a2ae208SSatish Balay #undef __FUNCT__ 16854a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 16869b94acceSBarry Smith /*@ 16879b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 168863a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 16899b94acceSBarry Smith 1690c7afd0dbSLois Curfman McInnes Collective on SNES 1691c7afd0dbSLois Curfman McInnes 1692b2002411SLois Curfman McInnes Input Parameters: 1693c7afd0dbSLois Curfman McInnes + snes - the SNES context 1694c7afd0dbSLois Curfman McInnes - x - the solution vector 16959b94acceSBarry Smith 1696b2002411SLois Curfman McInnes Notes: 16978ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 16988ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 16998ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 17008ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 17018ddd3da0SLois Curfman McInnes 170236851e7fSLois Curfman McInnes Level: beginner 170336851e7fSLois Curfman McInnes 17049b94acceSBarry Smith .keywords: SNES, nonlinear, solve 17059b94acceSBarry Smith 170663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 17079b94acceSBarry Smith @*/ 1708c9780b6fSBarry Smith int SNESSolve(SNES snes,Vec x) 17099b94acceSBarry Smith { 1710f1af5d2fSBarry Smith int ierr; 1711f1af5d2fSBarry Smith PetscTruth flg; 1712052efed2SBarry Smith 17133a40ed3dSBarry Smith PetscFunctionBegin; 17144482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17154482741eSBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1716c9780b6fSBarry Smith PetscCheckSameComm(snes,1,x,2); 171729bbc08cSBarry Smith if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 171874637425SBarry Smith 171982bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1720c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 1721758f92a0SBarry Smith if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 1722d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 172350ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1724c9780b6fSBarry Smith ierr = (*(snes)->solve)(snes);CHKERRQ(ierr); 1725d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1726b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 17278b179ff8SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1728da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1729da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 1730*5968eb51SBarry Smith if (snes->printreason) { 1731*5968eb51SBarry Smith if (snes->reason > 0) { 1732*5968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 1733*5968eb51SBarry Smith } else { 1734*5968eb51SBarry Smith ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 1735*5968eb51SBarry Smith } 1736*5968eb51SBarry Smith } 1737*5968eb51SBarry Smith 17383a40ed3dSBarry Smith PetscFunctionReturn(0); 17399b94acceSBarry Smith } 17409b94acceSBarry Smith 17419b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 17429b94acceSBarry Smith 17434a2ae208SSatish Balay #undef __FUNCT__ 17444a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 174582bf6240SBarry Smith /*@C 17464b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 17479b94acceSBarry Smith 1748fee21e36SBarry Smith Collective on SNES 1749fee21e36SBarry Smith 1750c7afd0dbSLois Curfman McInnes Input Parameters: 1751c7afd0dbSLois Curfman McInnes + snes - the SNES context 1752454a90a3SBarry Smith - type - a known method 1753c7afd0dbSLois Curfman McInnes 1754c7afd0dbSLois Curfman McInnes Options Database Key: 1755454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1756c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1757ae12b187SLois Curfman McInnes 17589b94acceSBarry Smith Notes: 1759e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 17604b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1761c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17624b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1763c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17649b94acceSBarry Smith 1765ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1766ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1767ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1768ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1769ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1770ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1771ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1772ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1773ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1774b0a32e0cSBarry Smith appropriate method. 177536851e7fSLois Curfman McInnes 177636851e7fSLois Curfman McInnes Level: intermediate 1777a703fe33SLois Curfman McInnes 1778454a90a3SBarry Smith .keywords: SNES, set, type 1779435da068SBarry Smith 1780435da068SBarry Smith .seealso: SNESType, SNESCreate() 1781435da068SBarry Smith 17829b94acceSBarry Smith @*/ 17830e33f6ddSBarry Smith int SNESSetType(SNES snes,const SNESType type) 17849b94acceSBarry Smith { 17850f5bd95cSBarry Smith int ierr,(*r)(SNES); 17866831982aSBarry Smith PetscTruth match; 17873a40ed3dSBarry Smith 17883a40ed3dSBarry Smith PetscFunctionBegin; 17894482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 17904482741eSBarry Smith PetscValidCharPointer(type,2); 179182bf6240SBarry Smith 17926831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 17930f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 179482bf6240SBarry Smith 179582bf6240SBarry Smith if (snes->setupcalled) { 1796e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 179782bf6240SBarry Smith snes->data = 0; 179882bf6240SBarry Smith } 17997d1a2b2bSBarry Smith 18009b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 180182bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 180282bf6240SBarry Smith 1803b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 180482bf6240SBarry Smith 180529bbc08cSBarry Smith if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type); 18061548b14aSBarry Smith 1807606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 180882bf6240SBarry Smith snes->data = 0; 18093a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1810454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 181182bf6240SBarry Smith 18123a40ed3dSBarry Smith PetscFunctionReturn(0); 18139b94acceSBarry Smith } 18149b94acceSBarry Smith 1815a847f771SSatish Balay 18169b94acceSBarry Smith /* --------------------------------------------------------------------- */ 18174a2ae208SSatish Balay #undef __FUNCT__ 18184a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 18199b94acceSBarry Smith /*@C 18209b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1821f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 18229b94acceSBarry Smith 1823fee21e36SBarry Smith Not Collective 1824fee21e36SBarry Smith 182536851e7fSLois Curfman McInnes Level: advanced 182636851e7fSLois Curfman McInnes 18279b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 18289b94acceSBarry Smith 18299b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 18309b94acceSBarry Smith @*/ 1831cf256101SBarry Smith int SNESRegisterDestroy(void) 18329b94acceSBarry Smith { 183382bf6240SBarry Smith int ierr; 183482bf6240SBarry Smith 18353a40ed3dSBarry Smith PetscFunctionBegin; 183682bf6240SBarry Smith if (SNESList) { 1837b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 183882bf6240SBarry Smith SNESList = 0; 18399b94acceSBarry Smith } 18404c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 18413a40ed3dSBarry Smith PetscFunctionReturn(0); 18429b94acceSBarry Smith } 18439b94acceSBarry Smith 18444a2ae208SSatish Balay #undef __FUNCT__ 18454a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 18469b94acceSBarry Smith /*@C 18479a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 18489b94acceSBarry Smith 1849c7afd0dbSLois Curfman McInnes Not Collective 1850c7afd0dbSLois Curfman McInnes 18519b94acceSBarry Smith Input Parameter: 18524b0e389bSBarry Smith . snes - nonlinear solver context 18539b94acceSBarry Smith 18549b94acceSBarry Smith Output Parameter: 18553a7fca6bSBarry Smith . type - SNES method (a character string) 18569b94acceSBarry Smith 185736851e7fSLois Curfman McInnes Level: intermediate 185836851e7fSLois Curfman McInnes 1859454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 18609b94acceSBarry Smith @*/ 1861454a90a3SBarry Smith int SNESGetType(SNES snes,SNESType *type) 18629b94acceSBarry Smith { 18633a40ed3dSBarry Smith PetscFunctionBegin; 18644482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18654482741eSBarry Smith PetscValidPointer(type,2); 1866454a90a3SBarry Smith *type = snes->type_name; 18673a40ed3dSBarry Smith PetscFunctionReturn(0); 18689b94acceSBarry Smith } 18699b94acceSBarry Smith 18704a2ae208SSatish Balay #undef __FUNCT__ 18714a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 18729b94acceSBarry Smith /*@C 18739b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 18749b94acceSBarry Smith stored. 18759b94acceSBarry Smith 1876c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1877c7afd0dbSLois Curfman McInnes 18789b94acceSBarry Smith Input Parameter: 18799b94acceSBarry Smith . snes - the SNES context 18809b94acceSBarry Smith 18819b94acceSBarry Smith Output Parameter: 18829b94acceSBarry Smith . x - the solution 18839b94acceSBarry Smith 188436851e7fSLois Curfman McInnes Level: advanced 188536851e7fSLois Curfman McInnes 18869b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 18879b94acceSBarry Smith 18884b27c08aSLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetSolutionUpdate() 18899b94acceSBarry Smith @*/ 18909b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 18919b94acceSBarry Smith { 18923a40ed3dSBarry Smith PetscFunctionBegin; 18934482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 18944482741eSBarry Smith PetscValidPointer(x,2); 18959b94acceSBarry Smith *x = snes->vec_sol_always; 18963a40ed3dSBarry Smith PetscFunctionReturn(0); 18979b94acceSBarry Smith } 18989b94acceSBarry Smith 18994a2ae208SSatish Balay #undef __FUNCT__ 19004a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 19019b94acceSBarry Smith /*@C 19029b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 19039b94acceSBarry Smith stored. 19049b94acceSBarry Smith 1905c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1906c7afd0dbSLois Curfman McInnes 19079b94acceSBarry Smith Input Parameter: 19089b94acceSBarry Smith . snes - the SNES context 19099b94acceSBarry Smith 19109b94acceSBarry Smith Output Parameter: 19119b94acceSBarry Smith . x - the solution update 19129b94acceSBarry Smith 191336851e7fSLois Curfman McInnes Level: advanced 191436851e7fSLois Curfman McInnes 19159b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 19169b94acceSBarry Smith 19179b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 19189b94acceSBarry Smith @*/ 19199b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 19209b94acceSBarry Smith { 19213a40ed3dSBarry Smith PetscFunctionBegin; 19224482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19234482741eSBarry Smith PetscValidPointer(x,2); 19249b94acceSBarry Smith *x = snes->vec_sol_update_always; 19253a40ed3dSBarry Smith PetscFunctionReturn(0); 19269b94acceSBarry Smith } 19279b94acceSBarry Smith 19284a2ae208SSatish Balay #undef __FUNCT__ 19294a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 19309b94acceSBarry Smith /*@C 19313638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 19329b94acceSBarry Smith 1933c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1934c7afd0dbSLois Curfman McInnes 19359b94acceSBarry Smith Input Parameter: 19369b94acceSBarry Smith . snes - the SNES context 19379b94acceSBarry Smith 19389b94acceSBarry Smith Output Parameter: 19397bf4e008SBarry Smith + r - the function (or PETSC_NULL) 194000036973SBarry Smith . ctx - the function context (or PETSC_NULL) 194100036973SBarry Smith - func - the function (or PETSC_NULL) 19429b94acceSBarry Smith 194336851e7fSLois Curfman McInnes Level: advanced 194436851e7fSLois Curfman McInnes 1945a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 19469b94acceSBarry Smith 19474b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 19489b94acceSBarry Smith @*/ 194900036973SBarry Smith int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*)) 19509b94acceSBarry Smith { 19513a40ed3dSBarry Smith PetscFunctionBegin; 19524482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 19537bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 19547bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 195500036973SBarry Smith if (func) *func = snes->computefunction; 19563a40ed3dSBarry Smith PetscFunctionReturn(0); 19579b94acceSBarry Smith } 19589b94acceSBarry Smith 19594a2ae208SSatish Balay #undef __FUNCT__ 19604a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 19613c7409f5SSatish Balay /*@C 19623c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1963d850072dSLois Curfman McInnes SNES options in the database. 19643c7409f5SSatish Balay 1965fee21e36SBarry Smith Collective on SNES 1966fee21e36SBarry Smith 1967c7afd0dbSLois Curfman McInnes Input Parameter: 1968c7afd0dbSLois Curfman McInnes + snes - the SNES context 1969c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1970c7afd0dbSLois Curfman McInnes 1971d850072dSLois Curfman McInnes Notes: 1972a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1973c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1974d850072dSLois Curfman McInnes 197536851e7fSLois Curfman McInnes Level: advanced 197636851e7fSLois Curfman McInnes 19773c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1978a86d99e1SLois Curfman McInnes 1979a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19803c7409f5SSatish Balay @*/ 19811836bdbcSSatish Balay int SNESSetOptionsPrefix(SNES snes,const char prefix[]) 19823c7409f5SSatish Balay { 19833c7409f5SSatish Balay int ierr; 19843c7409f5SSatish Balay 19853a40ed3dSBarry Smith PetscFunctionBegin; 19864482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1987639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 198894b7f48cSBarry Smith ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 19893a40ed3dSBarry Smith PetscFunctionReturn(0); 19903c7409f5SSatish Balay } 19913c7409f5SSatish Balay 19924a2ae208SSatish Balay #undef __FUNCT__ 19934a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 19943c7409f5SSatish Balay /*@C 1995f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1996d850072dSLois Curfman McInnes SNES options in the database. 19973c7409f5SSatish Balay 1998fee21e36SBarry Smith Collective on SNES 1999fee21e36SBarry Smith 2000c7afd0dbSLois Curfman McInnes Input Parameters: 2001c7afd0dbSLois Curfman McInnes + snes - the SNES context 2002c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2003c7afd0dbSLois Curfman McInnes 2004d850072dSLois Curfman McInnes Notes: 2005a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2006c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2007d850072dSLois Curfman McInnes 200836851e7fSLois Curfman McInnes Level: advanced 200936851e7fSLois Curfman McInnes 20103c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2011a86d99e1SLois Curfman McInnes 2012a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20133c7409f5SSatish Balay @*/ 20141836bdbcSSatish Balay int SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 20153c7409f5SSatish Balay { 20163c7409f5SSatish Balay int ierr; 20173c7409f5SSatish Balay 20183a40ed3dSBarry Smith PetscFunctionBegin; 20194482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2020639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 202194b7f48cSBarry Smith ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 20223a40ed3dSBarry Smith PetscFunctionReturn(0); 20233c7409f5SSatish Balay } 20243c7409f5SSatish Balay 20254a2ae208SSatish Balay #undef __FUNCT__ 20264a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 20279ab63eb5SSatish Balay /*@C 20283c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20293c7409f5SSatish Balay SNES options in the database. 20303c7409f5SSatish Balay 2031c7afd0dbSLois Curfman McInnes Not Collective 2032c7afd0dbSLois Curfman McInnes 20333c7409f5SSatish Balay Input Parameter: 20343c7409f5SSatish Balay . snes - the SNES context 20353c7409f5SSatish Balay 20363c7409f5SSatish Balay Output Parameter: 20373c7409f5SSatish Balay . prefix - pointer to the prefix string used 20383c7409f5SSatish Balay 20399ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 20409ab63eb5SSatish Balay sufficient length to hold the prefix. 20419ab63eb5SSatish Balay 204236851e7fSLois Curfman McInnes Level: advanced 204336851e7fSLois Curfman McInnes 20443c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2045a86d99e1SLois Curfman McInnes 2046a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20473c7409f5SSatish Balay @*/ 20481836bdbcSSatish Balay int SNESGetOptionsPrefix(SNES snes,char *prefix[]) 20493c7409f5SSatish Balay { 20503c7409f5SSatish Balay int ierr; 20513c7409f5SSatish Balay 20523a40ed3dSBarry Smith PetscFunctionBegin; 20534482741eSBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2054639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 20553a40ed3dSBarry Smith PetscFunctionReturn(0); 20563c7409f5SSatish Balay } 20573c7409f5SSatish Balay 2058b2002411SLois Curfman McInnes 20594a2ae208SSatish Balay #undef __FUNCT__ 20604a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 20613cea93caSBarry Smith /*@C 20623cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 20633cea93caSBarry Smith 20647f6c08e0SMatthew Knepley Level: advanced 20653cea93caSBarry Smith @*/ 20661836bdbcSSatish Balay int SNESRegister(const char sname[],const char path[],const char name[],int (*function)(SNES)) 2067b2002411SLois Curfman McInnes { 2068b2002411SLois Curfman McInnes char fullname[256]; 2069b2002411SLois Curfman McInnes int ierr; 2070b2002411SLois Curfman McInnes 2071b2002411SLois Curfman McInnes PetscFunctionBegin; 2072b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2073c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2074b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2075b2002411SLois Curfman McInnes } 2076da9b6338SBarry Smith 2077da9b6338SBarry Smith #undef __FUNCT__ 2078da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 2079da9b6338SBarry Smith int SNESTestLocalMin(SNES snes) 2080da9b6338SBarry Smith { 2081da9b6338SBarry Smith int ierr,N,i,j; 2082da9b6338SBarry Smith Vec u,uh,fh; 2083da9b6338SBarry Smith PetscScalar value; 2084da9b6338SBarry Smith PetscReal norm; 2085da9b6338SBarry Smith 2086da9b6338SBarry Smith PetscFunctionBegin; 2087da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2088da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2089da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2090da9b6338SBarry Smith 2091da9b6338SBarry Smith /* currently only works for sequential */ 2092da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2093da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2094da9b6338SBarry Smith for (i=0; i<N; i++) { 2095da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 2096da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %d\n",i);CHKERRQ(ierr); 2097da9b6338SBarry Smith for (j=-10; j<11; j++) { 2098ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2099da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2100da9b6338SBarry Smith ierr = (*snes->computefunction)(snes,uh,fh,snes->funP);CHKERRQ(ierr); 2101da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2102da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %d %18.16e\n",j,norm);CHKERRQ(ierr); 2103da9b6338SBarry Smith value = -value; 2104da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2105da9b6338SBarry Smith } 2106da9b6338SBarry Smith } 2107da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2108da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2109da9b6338SBarry Smith PetscFunctionReturn(0); 2110da9b6338SBarry Smith } 2111