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 */ 98ba1e511SMatthew Knepley int SNES_COOKIE; 104b27c08aSLois Curfman McInnes int SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval; 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; 479b94acceSBarry Smith SLES sles; 48454a90a3SBarry Smith char *type; 496831982aSBarry Smith PetscTruth isascii,isstring; 509b94acceSBarry Smith 513a40ed3dSBarry Smith PetscFunctionBegin; 5274679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 53b0a32e0cSBarry Smith if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm); 54b0a32e0cSBarry Smith PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE); 556831982aSBarry Smith PetscCheckSameComm(snes,viewer); 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 } 9377ed5343SBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 94b0a32e0cSBarry Smith ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 959b94acceSBarry Smith ierr = SLESView(sles,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" 11076b2cf59SMatthew Knepley /*@ 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 /*@ 136682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES 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 16236851e7fSLois Curfman McInnes - -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 16382738288SBarry Smith 16482738288SBarry Smith Options Database for Eisenstat-Walker method: 1654b27c08aSLois Curfman McInnes + -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence 1664b27c08aSLois Curfman McInnes . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 16736851e7fSLois Curfman McInnes . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 16836851e7fSLois Curfman McInnes . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 16936851e7fSLois Curfman McInnes . -snes_ksp_ew_gamma <gamma> - Sets gamma 17036851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha <alpha> - Sets alpha 17136851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 17236851e7fSLois Curfman McInnes - -snes_ksp_ew_threshold <threshold> - Sets threshold 17382738288SBarry Smith 17411ca99fdSLois Curfman McInnes Notes: 17511ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 17611ca99fdSLois Curfman McInnes the users manual. 17783e2fdc7SBarry Smith 17836851e7fSLois Curfman McInnes Level: beginner 17936851e7fSLois Curfman McInnes 1809b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1819b94acceSBarry Smith 18269ed319cSSatish Balay .seealso: SNESSetOptionsPrefix() 1839b94acceSBarry Smith @*/ 1849b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1859b94acceSBarry Smith { 1869b94acceSBarry Smith SLES sles; 187186905e3SBarry Smith SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 188f1af5d2fSBarry Smith PetscTruth flg; 18976b2cf59SMatthew Knepley int ierr, i; 190186905e3SBarry Smith char *deft,type[256]; 1919b94acceSBarry Smith 1923a40ed3dSBarry Smith PetscFunctionBegin; 19377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 194ca161407SBarry Smith 195b0a32e0cSBarry Smith ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 196186905e3SBarry Smith if (snes->type_name) { 197186905e3SBarry Smith deft = snes->type_name; 198186905e3SBarry Smith } else { 1994b27c08aSLois Curfman McInnes deft = SNESLS; 200d64ed03dSBarry Smith } 2014bbc92c1SBarry Smith 202186905e3SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 203b0a32e0cSBarry Smith ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 204d64ed03dSBarry Smith if (flg) { 205186905e3SBarry Smith ierr = SNESSetType(snes,type);CHKERRQ(ierr); 206186905e3SBarry Smith } else if (!snes->type_name) { 207186905e3SBarry Smith ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 208d64ed03dSBarry Smith } 209909c8a9fSBarry Smith ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 21093c39befSBarry Smith 21187828ca2SBarry Smith ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 21287828ca2SBarry Smith ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->atol,&snes->atol,0);CHKERRQ(ierr); 213186905e3SBarry Smith 21487828ca2SBarry Smith ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 215b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 216b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 21750ffb88aSMatthew Knepley ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 218186905e3SBarry Smith 219b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr); 220186905e3SBarry Smith 221b0a32e0cSBarry Smith ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 22287828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 22387828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 22487828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 22587828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 22687828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 22787828ca2SBarry Smith ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 228186905e3SBarry Smith 229b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr); 23093c39befSBarry Smith if (flg) {snes->converged = 0;} 231b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr); 2325cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 233b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr); 234b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);} 2353a7fca6bSBarry Smith ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr); 2363a7fca6bSBarry Smith if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);} 237b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr); 238b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);} 239b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr); 240b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 241b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr); 2427c922b88SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);} 2435ed2d596SBarry Smith ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr); 2445ed2d596SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);} 245b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr); 246186905e3SBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 247e24b481bSBarry Smith 248b0a32e0cSBarry Smith ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 2494b27c08aSLois Curfman McInnes if (flg) { 250186905e3SBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 251b0a32e0cSBarry Smith PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 2529b94acceSBarry Smith } 253639f9d9dSBarry Smith 25476b2cf59SMatthew Knepley for(i = 0; i < numberofsetfromoptions; i++) { 25576b2cf59SMatthew Knepley ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 25676b2cf59SMatthew Knepley } 25776b2cf59SMatthew Knepley 258186905e3SBarry Smith if (snes->setfromoptions) { 259186905e3SBarry Smith ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 260639f9d9dSBarry Smith } 261639f9d9dSBarry Smith 262b0a32e0cSBarry Smith ierr = PetscOptionsEnd();CHKERRQ(ierr); 2634bbc92c1SBarry Smith 2649b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 2659b94acceSBarry Smith ierr = SLESSetFromOptions(sles);CHKERRQ(ierr); 26693993e2dSLois Curfman McInnes 2673a40ed3dSBarry Smith PetscFunctionReturn(0); 2689b94acceSBarry Smith } 2699b94acceSBarry Smith 270a847f771SSatish Balay 2714a2ae208SSatish Balay #undef __FUNCT__ 2724a2ae208SSatish Balay #define __FUNCT__ "SNESSetApplicationContext" 2739b94acceSBarry Smith /*@ 2749b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 2759b94acceSBarry Smith the nonlinear solvers. 2769b94acceSBarry Smith 277fee21e36SBarry Smith Collective on SNES 278fee21e36SBarry Smith 279c7afd0dbSLois Curfman McInnes Input Parameters: 280c7afd0dbSLois Curfman McInnes + snes - the SNES context 281c7afd0dbSLois Curfman McInnes - usrP - optional user context 282c7afd0dbSLois Curfman McInnes 28336851e7fSLois Curfman McInnes Level: intermediate 28436851e7fSLois Curfman McInnes 2859b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 2869b94acceSBarry Smith 2879b94acceSBarry Smith .seealso: SNESGetApplicationContext() 2889b94acceSBarry Smith @*/ 2899b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 2909b94acceSBarry Smith { 2913a40ed3dSBarry Smith PetscFunctionBegin; 29277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2939b94acceSBarry Smith snes->user = usrP; 2943a40ed3dSBarry Smith PetscFunctionReturn(0); 2959b94acceSBarry Smith } 29674679c65SBarry Smith 2974a2ae208SSatish Balay #undef __FUNCT__ 2984a2ae208SSatish Balay #define __FUNCT__ "SNESGetApplicationContext" 2999b94acceSBarry Smith /*@C 3009b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3019b94acceSBarry Smith nonlinear solvers. 3029b94acceSBarry Smith 303c7afd0dbSLois Curfman McInnes Not Collective 304c7afd0dbSLois Curfman McInnes 3059b94acceSBarry Smith Input Parameter: 3069b94acceSBarry Smith . snes - SNES context 3079b94acceSBarry Smith 3089b94acceSBarry Smith Output Parameter: 3099b94acceSBarry Smith . usrP - user context 3109b94acceSBarry Smith 31136851e7fSLois Curfman McInnes Level: intermediate 31236851e7fSLois Curfman McInnes 3139b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3149b94acceSBarry Smith 3159b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3169b94acceSBarry Smith @*/ 3179b94acceSBarry Smith int SNESGetApplicationContext(SNES snes,void **usrP) 3189b94acceSBarry Smith { 3193a40ed3dSBarry Smith PetscFunctionBegin; 32077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3219b94acceSBarry Smith *usrP = snes->user; 3223a40ed3dSBarry Smith PetscFunctionReturn(0); 3239b94acceSBarry Smith } 32474679c65SBarry Smith 3254a2ae208SSatish Balay #undef __FUNCT__ 3264a2ae208SSatish Balay #define __FUNCT__ "SNESGetIterationNumber" 3279b94acceSBarry Smith /*@ 328c8228a4eSBarry Smith SNESGetIterationNumber - Gets the number of nonlinear iterations completed 329c8228a4eSBarry Smith at this time. 3309b94acceSBarry Smith 331c7afd0dbSLois Curfman McInnes Not Collective 332c7afd0dbSLois Curfman McInnes 3339b94acceSBarry Smith Input Parameter: 3349b94acceSBarry Smith . snes - SNES context 3359b94acceSBarry Smith 3369b94acceSBarry Smith Output Parameter: 3379b94acceSBarry Smith . iter - iteration number 3389b94acceSBarry Smith 339c8228a4eSBarry Smith Notes: 340c8228a4eSBarry Smith For example, during the computation of iteration 2 this would return 1. 341c8228a4eSBarry Smith 342c8228a4eSBarry Smith This is useful for using lagged Jacobians (where one does not recompute the 34308405cd6SLois Curfman McInnes Jacobian at each SNES iteration). For example, the code 34408405cd6SLois Curfman McInnes .vb 34508405cd6SLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&it); 34608405cd6SLois Curfman McInnes if (!(it % 2)) { 34708405cd6SLois Curfman McInnes [compute Jacobian here] 34808405cd6SLois Curfman McInnes } 34908405cd6SLois Curfman McInnes .ve 350c8228a4eSBarry Smith can be used in your ComputeJacobian() function to cause the Jacobian to be 35108405cd6SLois Curfman McInnes recomputed every second SNES iteration. 352c8228a4eSBarry Smith 35336851e7fSLois Curfman McInnes Level: intermediate 35436851e7fSLois Curfman McInnes 3559b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3569b94acceSBarry Smith @*/ 3579b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3589b94acceSBarry Smith { 3593a40ed3dSBarry Smith PetscFunctionBegin; 36077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 36174679c65SBarry Smith PetscValidIntPointer(iter); 3629b94acceSBarry Smith *iter = snes->iter; 3633a40ed3dSBarry Smith PetscFunctionReturn(0); 3649b94acceSBarry Smith } 36574679c65SBarry Smith 3664a2ae208SSatish Balay #undef __FUNCT__ 3674a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunctionNorm" 3689b94acceSBarry Smith /*@ 3699b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3709b94acceSBarry Smith with SNESSSetFunction(). 3719b94acceSBarry Smith 372c7afd0dbSLois Curfman McInnes Collective on SNES 373c7afd0dbSLois Curfman McInnes 3749b94acceSBarry Smith Input Parameter: 3759b94acceSBarry Smith . snes - SNES context 3769b94acceSBarry Smith 3779b94acceSBarry Smith Output Parameter: 3789b94acceSBarry Smith . fnorm - 2-norm of function 3799b94acceSBarry Smith 38036851e7fSLois Curfman McInnes Level: intermediate 38136851e7fSLois Curfman McInnes 3829b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 383a86d99e1SLois Curfman McInnes 38408405cd6SLois Curfman McInnes .seealso: SNESGetFunction() 3859b94acceSBarry Smith @*/ 38687828ca2SBarry Smith int SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm) 3879b94acceSBarry Smith { 3883a40ed3dSBarry Smith PetscFunctionBegin; 38977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 39074679c65SBarry Smith PetscValidScalarPointer(fnorm); 3919b94acceSBarry Smith *fnorm = snes->norm; 3923a40ed3dSBarry Smith PetscFunctionReturn(0); 3939b94acceSBarry Smith } 39474679c65SBarry Smith 3954a2ae208SSatish Balay #undef __FUNCT__ 3964a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps" 3979b94acceSBarry Smith /*@ 3989b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 3999b94acceSBarry Smith attempted by the nonlinear solver. 4009b94acceSBarry Smith 401c7afd0dbSLois Curfman McInnes Not Collective 402c7afd0dbSLois Curfman McInnes 4039b94acceSBarry Smith Input Parameter: 4049b94acceSBarry Smith . snes - SNES context 4059b94acceSBarry Smith 4069b94acceSBarry Smith Output Parameter: 4079b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4089b94acceSBarry Smith 409c96a6f78SLois Curfman McInnes Notes: 410c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 411c96a6f78SLois Curfman McInnes 41236851e7fSLois Curfman McInnes Level: intermediate 41336851e7fSLois Curfman McInnes 4149b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4159b94acceSBarry Smith @*/ 4169b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4179b94acceSBarry Smith { 4183a40ed3dSBarry Smith PetscFunctionBegin; 41977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 42074679c65SBarry Smith PetscValidIntPointer(nfails); 42150ffb88aSMatthew Knepley *nfails = snes->numFailures; 42250ffb88aSMatthew Knepley PetscFunctionReturn(0); 42350ffb88aSMatthew Knepley } 42450ffb88aSMatthew Knepley 42550ffb88aSMatthew Knepley #undef __FUNCT__ 42650ffb88aSMatthew Knepley #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps" 42750ffb88aSMatthew Knepley /*@ 42850ffb88aSMatthew Knepley SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps 42950ffb88aSMatthew Knepley attempted by the nonlinear solver before it gives up. 43050ffb88aSMatthew Knepley 43150ffb88aSMatthew Knepley Not Collective 43250ffb88aSMatthew Knepley 43350ffb88aSMatthew Knepley Input Parameters: 43450ffb88aSMatthew Knepley + snes - SNES context 43550ffb88aSMatthew Knepley - maxFails - maximum of unsuccessful steps 43650ffb88aSMatthew Knepley 43750ffb88aSMatthew Knepley Level: intermediate 43850ffb88aSMatthew Knepley 43950ffb88aSMatthew Knepley .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 44050ffb88aSMatthew Knepley @*/ 44150ffb88aSMatthew Knepley int SNESSetMaximumUnsuccessfulSteps(SNES snes, int maxFails) 44250ffb88aSMatthew Knepley { 44350ffb88aSMatthew Knepley PetscFunctionBegin; 44450ffb88aSMatthew Knepley PetscValidHeaderSpecific(snes,SNES_COOKIE); 44550ffb88aSMatthew Knepley snes->maxFailures = maxFails; 44650ffb88aSMatthew Knepley PetscFunctionReturn(0); 44750ffb88aSMatthew Knepley } 44850ffb88aSMatthew Knepley 44950ffb88aSMatthew Knepley #undef __FUNCT__ 45050ffb88aSMatthew Knepley #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps" 45150ffb88aSMatthew Knepley /*@ 45250ffb88aSMatthew Knepley SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps 45350ffb88aSMatthew Knepley attempted by the nonlinear solver before it gives up. 45450ffb88aSMatthew Knepley 45550ffb88aSMatthew Knepley Not Collective 45650ffb88aSMatthew Knepley 45750ffb88aSMatthew Knepley Input Parameter: 45850ffb88aSMatthew Knepley . snes - SNES context 45950ffb88aSMatthew Knepley 46050ffb88aSMatthew Knepley Output Parameter: 46150ffb88aSMatthew Knepley . maxFails - maximum of unsuccessful steps 46250ffb88aSMatthew Knepley 46350ffb88aSMatthew Knepley Level: intermediate 46450ffb88aSMatthew Knepley 46550ffb88aSMatthew Knepley .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 46650ffb88aSMatthew Knepley @*/ 46750ffb88aSMatthew Knepley int SNESGetMaximumUnsuccessfulSteps(SNES snes, int *maxFails) 46850ffb88aSMatthew Knepley { 46950ffb88aSMatthew Knepley PetscFunctionBegin; 47050ffb88aSMatthew Knepley PetscValidHeaderSpecific(snes,SNES_COOKIE); 47150ffb88aSMatthew Knepley PetscValidIntPointer(maxFails); 47250ffb88aSMatthew Knepley *maxFails = snes->maxFailures; 4733a40ed3dSBarry Smith PetscFunctionReturn(0); 4749b94acceSBarry Smith } 475a847f771SSatish Balay 4764a2ae208SSatish Balay #undef __FUNCT__ 4774a2ae208SSatish Balay #define __FUNCT__ "SNESGetNumberLinearIterations" 478c96a6f78SLois Curfman McInnes /*@ 479c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 480c96a6f78SLois Curfman McInnes used by the nonlinear solver. 481c96a6f78SLois Curfman McInnes 482c7afd0dbSLois Curfman McInnes Not Collective 483c7afd0dbSLois Curfman McInnes 484c96a6f78SLois Curfman McInnes Input Parameter: 485c96a6f78SLois Curfman McInnes . snes - SNES context 486c96a6f78SLois Curfman McInnes 487c96a6f78SLois Curfman McInnes Output Parameter: 488c96a6f78SLois Curfman McInnes . lits - number of linear iterations 489c96a6f78SLois Curfman McInnes 490c96a6f78SLois Curfman McInnes Notes: 491c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 492c96a6f78SLois Curfman McInnes 49336851e7fSLois Curfman McInnes Level: intermediate 49436851e7fSLois Curfman McInnes 495c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 496c96a6f78SLois Curfman McInnes @*/ 49786bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 498c96a6f78SLois Curfman McInnes { 4993a40ed3dSBarry Smith PetscFunctionBegin; 500c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 501c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 502c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 5033a40ed3dSBarry Smith PetscFunctionReturn(0); 504c96a6f78SLois Curfman McInnes } 505c96a6f78SLois Curfman McInnes 5064a2ae208SSatish Balay #undef __FUNCT__ 5074a2ae208SSatish Balay #define __FUNCT__ "SNESGetSLES" 5089b94acceSBarry Smith /*@C 5099b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 5109b94acceSBarry Smith 511c7afd0dbSLois Curfman McInnes Not Collective, but if SNES object is parallel, then SLES object is parallel 512c7afd0dbSLois Curfman McInnes 5139b94acceSBarry Smith Input Parameter: 5149b94acceSBarry Smith . snes - the SNES context 5159b94acceSBarry Smith 5169b94acceSBarry Smith Output Parameter: 5179b94acceSBarry Smith . sles - the SLES context 5189b94acceSBarry Smith 5199b94acceSBarry Smith Notes: 5209b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 5219b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5229b94acceSBarry Smith KSP and PC contexts as well. 5239b94acceSBarry Smith 52436851e7fSLois Curfman McInnes Level: beginner 52536851e7fSLois Curfman McInnes 5269b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 5279b94acceSBarry Smith 5289b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 5299b94acceSBarry Smith @*/ 5309b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 5319b94acceSBarry Smith { 5323a40ed3dSBarry Smith PetscFunctionBegin; 53377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5349b94acceSBarry Smith *sles = snes->sles; 5353a40ed3dSBarry Smith PetscFunctionReturn(0); 5369b94acceSBarry Smith } 53782bf6240SBarry Smith 5384a2ae208SSatish Balay #undef __FUNCT__ 5394a2ae208SSatish Balay #define __FUNCT__ "SNESPublish_Petsc" 540454a90a3SBarry Smith static int SNESPublish_Petsc(PetscObject obj) 541e24b481bSBarry Smith { 542aa482453SBarry Smith #if defined(PETSC_HAVE_AMS) 543454a90a3SBarry Smith SNES v = (SNES) obj; 544e24b481bSBarry Smith int ierr; 54543d6d2cbSBarry Smith #endif 546e24b481bSBarry Smith 547e24b481bSBarry Smith PetscFunctionBegin; 548e24b481bSBarry Smith 54943d6d2cbSBarry Smith #if defined(PETSC_HAVE_AMS) 550e24b481bSBarry Smith /* if it is already published then return */ 551e24b481bSBarry Smith if (v->amem >=0) PetscFunctionReturn(0); 552e24b481bSBarry Smith 553454a90a3SBarry Smith ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr); 554e24b481bSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ, 555e24b481bSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 556e24b481bSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ, 557e24b481bSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 558454a90a3SBarry Smith ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr); 559433994e6SBarry Smith #endif 560e24b481bSBarry Smith PetscFunctionReturn(0); 561e24b481bSBarry Smith } 562e24b481bSBarry Smith 5639b94acceSBarry Smith /* -----------------------------------------------------------*/ 5644a2ae208SSatish Balay #undef __FUNCT__ 5654a2ae208SSatish Balay #define __FUNCT__ "SNESCreate" 5669b94acceSBarry Smith /*@C 5679b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5689b94acceSBarry Smith 569c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 570c7afd0dbSLois Curfman McInnes 571c7afd0dbSLois Curfman McInnes Input Parameters: 572c7afd0dbSLois Curfman McInnes + comm - MPI communicator 5739b94acceSBarry Smith 5749b94acceSBarry Smith Output Parameter: 5759b94acceSBarry Smith . outsnes - the new SNES context 5769b94acceSBarry Smith 577c7afd0dbSLois Curfman McInnes Options Database Keys: 578c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 579c7afd0dbSLois Curfman McInnes and no preconditioning matrix 580c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 581c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 582c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 583c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 584c1f60f51SBarry Smith 58536851e7fSLois Curfman McInnes Level: beginner 58636851e7fSLois Curfman McInnes 5879b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5889b94acceSBarry Smith 5894b27c08aSLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy(), SNES 5909b94acceSBarry Smith @*/ 5914b27c08aSLois Curfman McInnes int SNESCreate(MPI_Comm comm,SNES *outsnes) 5929b94acceSBarry Smith { 5939b94acceSBarry Smith int ierr; 5949b94acceSBarry Smith SNES snes; 5959b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 59637fcc0dbSBarry Smith 5973a40ed3dSBarry Smith PetscFunctionBegin; 5988ba1e511SMatthew Knepley PetscValidPointer(outsnes); 5998ba1e511SMatthew Knepley *outsnes = PETSC_NULL; 6008ba1e511SMatthew Knepley #ifndef PETSC_USE_DYNAMIC_LIBRARIES 6018ba1e511SMatthew Knepley ierr = SNESInitializePackage(PETSC_NULL); CHKERRQ(ierr); 6028ba1e511SMatthew Knepley #endif 6038ba1e511SMatthew Knepley 6043f1db9ecSBarry Smith PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView); 605b0a32e0cSBarry Smith PetscLogObjectCreate(snes); 606e24b481bSBarry Smith snes->bops->publish = SNESPublish_Petsc; 6079b94acceSBarry Smith snes->max_its = 50; 6089750a799SBarry Smith snes->max_funcs = 10000; 6099b94acceSBarry Smith snes->norm = 0.0; 610b4874afaSBarry Smith snes->rtol = 1.e-8; 611b4874afaSBarry Smith snes->ttol = 0.0; 612b4874afaSBarry Smith snes->atol = 1.e-50; 6139b94acceSBarry Smith snes->xtol = 1.e-8; 6144b27c08aSLois Curfman McInnes snes->deltatol = 1.e-12; 6159b94acceSBarry Smith snes->nfuncs = 0; 61650ffb88aSMatthew Knepley snes->numFailures = 0; 61750ffb88aSMatthew Knepley snes->maxFailures = 1; 6187a00f4a9SLois Curfman McInnes snes->linear_its = 0; 619639f9d9dSBarry Smith snes->numbermonitors = 0; 6209b94acceSBarry Smith snes->data = 0; 6219b94acceSBarry Smith snes->view = 0; 62282bf6240SBarry Smith snes->setupcalled = 0; 623186905e3SBarry Smith snes->ksp_ewconv = PETSC_FALSE; 6246f24a144SLois Curfman McInnes snes->vwork = 0; 6256f24a144SLois Curfman McInnes snes->nwork = 0; 626758f92a0SBarry Smith snes->conv_hist_len = 0; 627758f92a0SBarry Smith snes->conv_hist_max = 0; 628758f92a0SBarry Smith snes->conv_hist = PETSC_NULL; 629758f92a0SBarry Smith snes->conv_hist_its = PETSC_NULL; 630758f92a0SBarry Smith snes->conv_hist_reset = PETSC_TRUE; 631184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 6329b94acceSBarry Smith 6339b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 634b0a32e0cSBarry Smith ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 635b0a32e0cSBarry Smith PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 6369b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6379b94acceSBarry Smith kctx->version = 2; 6389b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6399b94acceSBarry Smith this was too large for some test cases */ 6409b94acceSBarry Smith kctx->rtol_last = 0; 6419b94acceSBarry Smith kctx->rtol_max = .9; 6429b94acceSBarry Smith kctx->gamma = 1.0; 6439b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6449b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6459b94acceSBarry Smith kctx->threshold = .1; 6469b94acceSBarry Smith kctx->lresid_last = 0; 6479b94acceSBarry Smith kctx->norm_last = 0; 6489b94acceSBarry Smith 6499b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles);CHKERRQ(ierr); 650b0a32e0cSBarry Smith PetscLogObjectParent(snes,snes->sles) 6515334005bSBarry Smith 6529b94acceSBarry Smith *outsnes = snes; 65300036973SBarry Smith ierr = PetscPublishAll(snes);CHKERRQ(ierr); 6543a40ed3dSBarry Smith PetscFunctionReturn(0); 6559b94acceSBarry Smith } 6569b94acceSBarry Smith 6579b94acceSBarry Smith /* --------------------------------------------------------------- */ 6584a2ae208SSatish Balay #undef __FUNCT__ 6594a2ae208SSatish Balay #define __FUNCT__ "SNESSetFunction" 6609b94acceSBarry Smith /*@C 6619b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6629b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6639b94acceSBarry Smith equations. 6649b94acceSBarry Smith 665fee21e36SBarry Smith Collective on SNES 666fee21e36SBarry Smith 667c7afd0dbSLois Curfman McInnes Input Parameters: 668c7afd0dbSLois Curfman McInnes + snes - the SNES context 669c7afd0dbSLois Curfman McInnes . func - function evaluation routine 670c7afd0dbSLois Curfman McInnes . r - vector to store function value 671c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 672c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6739b94acceSBarry Smith 674c7afd0dbSLois Curfman McInnes Calling sequence of func: 6758d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 676c7afd0dbSLois Curfman McInnes 677313e4042SLois Curfman McInnes . f - function vector 678c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 6799b94acceSBarry Smith 6809b94acceSBarry Smith Notes: 6819b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6829b94acceSBarry Smith $ f'(x) x = -f(x), 683c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 6849b94acceSBarry Smith 68536851e7fSLois Curfman McInnes Level: beginner 68636851e7fSLois Curfman McInnes 6879b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6889b94acceSBarry Smith 689a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6909b94acceSBarry Smith @*/ 6915334005bSBarry Smith int SNESSetFunction(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 6929b94acceSBarry Smith { 6933a40ed3dSBarry Smith PetscFunctionBegin; 69477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 695184914b5SBarry Smith PetscValidHeaderSpecific(r,VEC_COOKIE); 696184914b5SBarry Smith PetscCheckSameComm(snes,r); 697184914b5SBarry Smith 6989b94acceSBarry Smith snes->computefunction = func; 6999b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7009b94acceSBarry Smith snes->funP = ctx; 7013a40ed3dSBarry Smith PetscFunctionReturn(0); 7029b94acceSBarry Smith } 7039b94acceSBarry Smith 7044a2ae208SSatish Balay #undef __FUNCT__ 7054a2ae208SSatish Balay #define __FUNCT__ "SNESComputeFunction" 7069b94acceSBarry Smith /*@ 70736851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7089b94acceSBarry Smith SNESSetFunction(). 7099b94acceSBarry Smith 710c7afd0dbSLois Curfman McInnes Collective on SNES 711c7afd0dbSLois Curfman McInnes 7129b94acceSBarry Smith Input Parameters: 713c7afd0dbSLois Curfman McInnes + snes - the SNES context 714c7afd0dbSLois Curfman McInnes - x - input vector 7159b94acceSBarry Smith 7169b94acceSBarry Smith Output Parameter: 7173638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7189b94acceSBarry Smith 7191bffabb2SLois Curfman McInnes Notes: 72036851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 72136851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 72236851e7fSLois Curfman McInnes themselves. 72336851e7fSLois Curfman McInnes 72436851e7fSLois Curfman McInnes Level: developer 72536851e7fSLois Curfman McInnes 7269b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 7279b94acceSBarry Smith 728a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 7299b94acceSBarry Smith @*/ 7309b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x,Vec y) 7319b94acceSBarry Smith { 7329b94acceSBarry Smith int ierr; 7339b94acceSBarry Smith 7343a40ed3dSBarry Smith PetscFunctionBegin; 735184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 736184914b5SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 737184914b5SBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE); 738184914b5SBarry Smith PetscCheckSameComm(snes,x); 739184914b5SBarry Smith PetscCheckSameComm(snes,y); 740184914b5SBarry Smith 741d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 742d64ed03dSBarry Smith PetscStackPush("SNES user function"); 7439b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 744d64ed03dSBarry Smith PetscStackPop; 745ae3c334cSLois Curfman McInnes snes->nfuncs++; 746d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 7473a40ed3dSBarry Smith PetscFunctionReturn(0); 7489b94acceSBarry Smith } 7499b94acceSBarry Smith 7504a2ae208SSatish Balay #undef __FUNCT__ 7514a2ae208SSatish Balay #define __FUNCT__ "SNESComputeJacobian" 75262fef451SLois Curfman McInnes /*@ 75362fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 75462fef451SLois Curfman McInnes set with SNESSetJacobian(). 75562fef451SLois Curfman McInnes 756c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 757c7afd0dbSLois Curfman McInnes 75862fef451SLois Curfman McInnes Input Parameters: 759c7afd0dbSLois Curfman McInnes + snes - the SNES context 760c7afd0dbSLois Curfman McInnes - x - input vector 76162fef451SLois Curfman McInnes 76262fef451SLois Curfman McInnes Output Parameters: 763c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 76462fef451SLois Curfman McInnes . B - optional preconditioning matrix 765c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 766fee21e36SBarry Smith 76762fef451SLois Curfman McInnes Notes: 76862fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 76962fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 77062fef451SLois Curfman McInnes 771dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 772dc5a77f8SLois Curfman McInnes flag parameter. 77362fef451SLois Curfman McInnes 77436851e7fSLois Curfman McInnes Level: developer 77536851e7fSLois Curfman McInnes 77662fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 77762fef451SLois Curfman McInnes 77862fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 77962fef451SLois Curfman McInnes @*/ 7809b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 7819b94acceSBarry Smith { 7829b94acceSBarry Smith int ierr; 7833a40ed3dSBarry Smith 7843a40ed3dSBarry Smith PetscFunctionBegin; 785184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 786184914b5SBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE); 787184914b5SBarry Smith PetscCheckSameComm(snes,X); 7883a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 789d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 790c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 791d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 7929b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 793d64ed03dSBarry Smith PetscStackPop; 794d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 7956d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 79677c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 79777c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 7983a40ed3dSBarry Smith PetscFunctionReturn(0); 7999b94acceSBarry Smith } 8009b94acceSBarry Smith 8014a2ae208SSatish Balay #undef __FUNCT__ 8024a2ae208SSatish Balay #define __FUNCT__ "SNESSetJacobian" 8039b94acceSBarry Smith /*@C 8049b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 805044dda88SLois Curfman McInnes location to store the matrix. 8069b94acceSBarry Smith 807c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 808c7afd0dbSLois Curfman McInnes 8099b94acceSBarry Smith Input Parameters: 810c7afd0dbSLois Curfman McInnes + snes - the SNES context 8119b94acceSBarry Smith . A - Jacobian matrix 8129b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 8139b94acceSBarry Smith . func - Jacobian evaluation routine 814c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 8152cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 8169b94acceSBarry Smith 8179b94acceSBarry Smith Calling sequence of func: 8188d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 8199b94acceSBarry Smith 820c7afd0dbSLois Curfman McInnes + x - input vector 8219b94acceSBarry Smith . A - Jacobian matrix 8229b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 823ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 824ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 825c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 8269b94acceSBarry Smith 8279b94acceSBarry Smith Notes: 828dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 8292cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 830ac21db08SLois Curfman McInnes 831ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 8329b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 8339b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 8349b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 8359b94acceSBarry Smith throughout the global iterations. 8369b94acceSBarry Smith 83736851e7fSLois Curfman McInnes Level: beginner 83836851e7fSLois Curfman McInnes 8399b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 8409b94acceSBarry Smith 841ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 8429b94acceSBarry Smith @*/ 843454a90a3SBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 8449b94acceSBarry Smith { 8453a7fca6bSBarry Smith int ierr; 8463a7fca6bSBarry Smith 8473a40ed3dSBarry Smith PetscFunctionBegin; 84877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 84900036973SBarry Smith if (A) PetscValidHeaderSpecific(A,MAT_COOKIE); 85000036973SBarry Smith if (B) PetscValidHeaderSpecific(B,MAT_COOKIE); 85100036973SBarry Smith if (A) PetscCheckSameComm(snes,A); 85200036973SBarry Smith if (B) PetscCheckSameComm(snes,B); 8533a7fca6bSBarry Smith if (func) snes->computejacobian = func; 8543a7fca6bSBarry Smith if (ctx) snes->jacP = ctx; 8553a7fca6bSBarry Smith if (A) { 8563a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 8579b94acceSBarry Smith snes->jacobian = A; 8583a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 8593a7fca6bSBarry Smith } 8603a7fca6bSBarry Smith if (B) { 8613a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 8629b94acceSBarry Smith snes->jacobian_pre = B; 8633a7fca6bSBarry Smith ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 8643a7fca6bSBarry Smith } 8653a40ed3dSBarry Smith PetscFunctionReturn(0); 8669b94acceSBarry Smith } 86762fef451SLois Curfman McInnes 8684a2ae208SSatish Balay #undef __FUNCT__ 8694a2ae208SSatish Balay #define __FUNCT__ "SNESGetJacobian" 870c2aafc4cSSatish Balay /*@C 871b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 872b4fd4287SBarry Smith provided context for evaluating the Jacobian. 873b4fd4287SBarry Smith 874c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 875c7afd0dbSLois Curfman McInnes 876b4fd4287SBarry Smith Input Parameter: 877b4fd4287SBarry Smith . snes - the nonlinear solver context 878b4fd4287SBarry Smith 879b4fd4287SBarry Smith Output Parameters: 880c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 881b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 88200036973SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 88300036973SBarry Smith - func - location to put Jacobian function (or PETSC_NULL) 884fee21e36SBarry Smith 88536851e7fSLois Curfman McInnes Level: advanced 88636851e7fSLois Curfman McInnes 887b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 888b4fd4287SBarry Smith @*/ 88900036973SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 890b4fd4287SBarry Smith { 8913a40ed3dSBarry Smith PetscFunctionBegin; 892184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 893b4fd4287SBarry Smith if (A) *A = snes->jacobian; 894b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 895b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 89600036973SBarry Smith if (func) *func = snes->computejacobian; 8973a40ed3dSBarry Smith PetscFunctionReturn(0); 898b4fd4287SBarry Smith } 899b4fd4287SBarry Smith 9009b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 90145fc7adcSBarry Smith extern int SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 9029b94acceSBarry Smith 9034a2ae208SSatish Balay #undef __FUNCT__ 9044a2ae208SSatish Balay #define __FUNCT__ "SNESSetUp" 9059b94acceSBarry Smith /*@ 9069b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 907272ac6f2SLois Curfman McInnes of a nonlinear solver. 9089b94acceSBarry Smith 909fee21e36SBarry Smith Collective on SNES 910fee21e36SBarry Smith 911c7afd0dbSLois Curfman McInnes Input Parameters: 912c7afd0dbSLois Curfman McInnes + snes - the SNES context 913c7afd0dbSLois Curfman McInnes - x - the solution vector 914c7afd0dbSLois Curfman McInnes 915272ac6f2SLois Curfman McInnes Notes: 916272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 917272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 918272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 919272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 920272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 921272ac6f2SLois Curfman McInnes 92236851e7fSLois Curfman McInnes Level: advanced 92336851e7fSLois Curfman McInnes 9249b94acceSBarry Smith .keywords: SNES, nonlinear, setup 9259b94acceSBarry Smith 9269b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 9279b94acceSBarry Smith @*/ 9288ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 9299b94acceSBarry Smith { 930f1af5d2fSBarry Smith int ierr; 9314b27c08aSLois Curfman McInnes PetscTruth flg, iseqtr; 9323a40ed3dSBarry Smith 9333a40ed3dSBarry Smith PetscFunctionBegin; 93477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 93577c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 936184914b5SBarry Smith PetscCheckSameComm(snes,x); 9378ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 9389b94acceSBarry Smith 939b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 940c1f60f51SBarry Smith /* 941c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 942dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 943c1f60f51SBarry Smith */ 944c1f60f51SBarry Smith if (flg) { 945c1f60f51SBarry Smith Mat J; 9465a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 9475a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 9483a7fca6bSBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n"); 9493a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 9503a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 951c1f60f51SBarry Smith } 95245fc7adcSBarry Smith 953b4014bb2SMatthew Knepley #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) 95445fc7adcSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 95545fc7adcSBarry Smith if (flg) { 95645fc7adcSBarry Smith Mat J; 95745fc7adcSBarry Smith ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 95845fc7adcSBarry Smith ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 95945fc7adcSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 96045fc7adcSBarry Smith } 96132a4b47aSMatthew Knepley #endif 96245fc7adcSBarry Smith 963b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 964c1f60f51SBarry Smith /* 965dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 966c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 967c1f60f51SBarry Smith */ 96831615d04SLois Curfman McInnes if (flg) { 969272ac6f2SLois Curfman McInnes Mat J; 970f3ef73ceSBarry Smith SLES sles; 971f3ef73ceSBarry Smith PC pc; 972f3ef73ceSBarry Smith 9735a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 9743a7fca6bSBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 97593b98531SBarry Smith PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"); 9763a7fca6bSBarry Smith ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 9773a7fca6bSBarry Smith ierr = MatDestroy(J);CHKERRQ(ierr); 9783a7fca6bSBarry Smith 979f3ef73ceSBarry Smith /* force no preconditioner */ 980f3ef73ceSBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 981f3ef73ceSBarry Smith ierr = SLESGetPC(sles,&pc);CHKERRQ(ierr); 982a9815358SBarry Smith ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 983a9815358SBarry Smith if (!flg) { 984f3ef73ceSBarry Smith ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 985272ac6f2SLois Curfman McInnes } 986a9815358SBarry Smith } 987f3ef73ceSBarry Smith 98829bbc08cSBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 98929bbc08cSBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 99029bbc08cSBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 991a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 99229bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 993a8c6a408SBarry Smith } 994a703fe33SLois Curfman McInnes 995a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 9964b27c08aSLois Curfman McInnes ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 9976831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 998a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 999a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 1000a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 1001186905e3SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1002a703fe33SLois Curfman McInnes } 10034b27c08aSLois Curfman McInnes 1004a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 100582bf6240SBarry Smith snes->setupcalled = 1; 10063a40ed3dSBarry Smith PetscFunctionReturn(0); 10079b94acceSBarry Smith } 10089b94acceSBarry Smith 10094a2ae208SSatish Balay #undef __FUNCT__ 10104a2ae208SSatish Balay #define __FUNCT__ "SNESDestroy" 10119b94acceSBarry Smith /*@C 10129b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 10139b94acceSBarry Smith with SNESCreate(). 10149b94acceSBarry Smith 1015c7afd0dbSLois Curfman McInnes Collective on SNES 1016c7afd0dbSLois Curfman McInnes 10179b94acceSBarry Smith Input Parameter: 10189b94acceSBarry Smith . snes - the SNES context 10199b94acceSBarry Smith 102036851e7fSLois Curfman McInnes Level: beginner 102136851e7fSLois Curfman McInnes 10229b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 10239b94acceSBarry Smith 102463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 10259b94acceSBarry Smith @*/ 10269b94acceSBarry Smith int SNESDestroy(SNES snes) 10279b94acceSBarry Smith { 1028b8a78c4aSBarry Smith int i,ierr; 10293a40ed3dSBarry Smith 10303a40ed3dSBarry Smith PetscFunctionBegin; 103177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 10323a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1033d4bb536fSBarry Smith 1034be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 10350f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1036be0abb6dSBarry Smith 1037e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1038606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 10393a7fca6bSBarry Smith if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 10403a7fca6bSBarry Smith if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 10419b94acceSBarry Smith ierr = SLESDestroy(snes->sles);CHKERRQ(ierr); 1042522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1043b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++) { 1044b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1045b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1046b8a78c4aSBarry Smith } 1047b8a78c4aSBarry Smith } 1048b0a32e0cSBarry Smith PetscLogObjectDestroy((PetscObject)snes); 10490452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 10503a40ed3dSBarry Smith PetscFunctionReturn(0); 10519b94acceSBarry Smith } 10529b94acceSBarry Smith 10539b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 10549b94acceSBarry Smith 10554a2ae208SSatish Balay #undef __FUNCT__ 10564a2ae208SSatish Balay #define __FUNCT__ "SNESSetTolerances" 10579b94acceSBarry Smith /*@ 1058d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 10599b94acceSBarry Smith 1060c7afd0dbSLois Curfman McInnes Collective on SNES 1061c7afd0dbSLois Curfman McInnes 10629b94acceSBarry Smith Input Parameters: 1063c7afd0dbSLois Curfman McInnes + snes - the SNES context 106433174efeSLois Curfman McInnes . atol - absolute convergence tolerance 106533174efeSLois Curfman McInnes . rtol - relative convergence tolerance 106633174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 106733174efeSLois Curfman McInnes of the change in the solution between steps 106833174efeSLois Curfman McInnes . maxit - maximum number of iterations 1069c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1070fee21e36SBarry Smith 107133174efeSLois Curfman McInnes Options Database Keys: 1072c7afd0dbSLois Curfman McInnes + -snes_atol <atol> - Sets atol 1073c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1074c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1075c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1076c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 10779b94acceSBarry Smith 1078d7a720efSLois Curfman McInnes Notes: 10799b94acceSBarry Smith The default maximum number of iterations is 50. 10809b94acceSBarry Smith The default maximum number of function evaluations is 1000. 10819b94acceSBarry Smith 108236851e7fSLois Curfman McInnes Level: intermediate 108336851e7fSLois Curfman McInnes 108433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 10859b94acceSBarry Smith 1086d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 10879b94acceSBarry Smith @*/ 1088329f5518SBarry Smith int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf) 10899b94acceSBarry Smith { 10903a40ed3dSBarry Smith PetscFunctionBegin; 109177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1092d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1093d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1094d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1095d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1096d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 10973a40ed3dSBarry Smith PetscFunctionReturn(0); 10989b94acceSBarry Smith } 10999b94acceSBarry Smith 11004a2ae208SSatish Balay #undef __FUNCT__ 11014a2ae208SSatish Balay #define __FUNCT__ "SNESGetTolerances" 11029b94acceSBarry Smith /*@ 110333174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 110433174efeSLois Curfman McInnes 1105c7afd0dbSLois Curfman McInnes Not Collective 1106c7afd0dbSLois Curfman McInnes 110733174efeSLois Curfman McInnes Input Parameters: 1108c7afd0dbSLois Curfman McInnes + snes - the SNES context 110933174efeSLois Curfman McInnes . atol - absolute convergence tolerance 111033174efeSLois Curfman McInnes . rtol - relative convergence tolerance 111133174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 111233174efeSLois Curfman McInnes of the change in the solution between steps 111333174efeSLois Curfman McInnes . maxit - maximum number of iterations 1114c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1115fee21e36SBarry Smith 111633174efeSLois Curfman McInnes Notes: 111733174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 111833174efeSLois Curfman McInnes 111936851e7fSLois Curfman McInnes Level: intermediate 112036851e7fSLois Curfman McInnes 112133174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 112233174efeSLois Curfman McInnes 112333174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 112433174efeSLois Curfman McInnes @*/ 1125329f5518SBarry Smith int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf) 112633174efeSLois Curfman McInnes { 11273a40ed3dSBarry Smith PetscFunctionBegin; 112833174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 112933174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 113033174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 113133174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 113233174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 113333174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 11343a40ed3dSBarry Smith PetscFunctionReturn(0); 113533174efeSLois Curfman McInnes } 113633174efeSLois Curfman McInnes 11374a2ae208SSatish Balay #undef __FUNCT__ 11384a2ae208SSatish Balay #define __FUNCT__ "SNESSetTrustRegionTolerance" 113933174efeSLois Curfman McInnes /*@ 11409b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 11419b94acceSBarry Smith 1142fee21e36SBarry Smith Collective on SNES 1143fee21e36SBarry Smith 1144c7afd0dbSLois Curfman McInnes Input Parameters: 1145c7afd0dbSLois Curfman McInnes + snes - the SNES context 1146c7afd0dbSLois Curfman McInnes - tol - tolerance 1147c7afd0dbSLois Curfman McInnes 11489b94acceSBarry Smith Options Database Key: 1149c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 11509b94acceSBarry Smith 115136851e7fSLois Curfman McInnes Level: intermediate 115236851e7fSLois Curfman McInnes 11539b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 11549b94acceSBarry Smith 1155d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 11569b94acceSBarry Smith @*/ 1157329f5518SBarry Smith int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 11589b94acceSBarry Smith { 11593a40ed3dSBarry Smith PetscFunctionBegin; 116077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 11619b94acceSBarry Smith snes->deltatol = tol; 11623a40ed3dSBarry Smith PetscFunctionReturn(0); 11639b94acceSBarry Smith } 11649b94acceSBarry Smith 1165df9fa365SBarry Smith /* 1166df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1167df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1168df9fa365SBarry Smith macros instead of functions 1169df9fa365SBarry Smith */ 11704a2ae208SSatish Balay #undef __FUNCT__ 11714a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitor" 1172329f5518SBarry Smith int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx) 1173ce1608b8SBarry Smith { 1174ce1608b8SBarry Smith int ierr; 1175ce1608b8SBarry Smith 1176ce1608b8SBarry Smith PetscFunctionBegin; 1177184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1178ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1179ce1608b8SBarry Smith PetscFunctionReturn(0); 1180ce1608b8SBarry Smith } 1181ce1608b8SBarry Smith 11824a2ae208SSatish Balay #undef __FUNCT__ 11834a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorCreate" 1184b0a32e0cSBarry Smith int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n,PetscDrawLG *draw) 1185df9fa365SBarry Smith { 1186df9fa365SBarry Smith int ierr; 1187df9fa365SBarry Smith 1188df9fa365SBarry Smith PetscFunctionBegin; 1189df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1190df9fa365SBarry Smith PetscFunctionReturn(0); 1191df9fa365SBarry Smith } 1192df9fa365SBarry Smith 11934a2ae208SSatish Balay #undef __FUNCT__ 11944a2ae208SSatish Balay #define __FUNCT__ "SNESLGMonitorDestroy" 1195b0a32e0cSBarry Smith int SNESLGMonitorDestroy(PetscDrawLG draw) 1196df9fa365SBarry Smith { 1197df9fa365SBarry Smith int ierr; 1198df9fa365SBarry Smith 1199df9fa365SBarry Smith PetscFunctionBegin; 1200df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1201df9fa365SBarry Smith PetscFunctionReturn(0); 1202df9fa365SBarry Smith } 1203df9fa365SBarry Smith 12049b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12059b94acceSBarry Smith 12064a2ae208SSatish Balay #undef __FUNCT__ 12074a2ae208SSatish Balay #define __FUNCT__ "SNESSetMonitor" 12089b94acceSBarry Smith /*@C 12095cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 12109b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12119b94acceSBarry Smith progress. 12129b94acceSBarry Smith 1213fee21e36SBarry Smith Collective on SNES 1214fee21e36SBarry Smith 1215c7afd0dbSLois Curfman McInnes Input Parameters: 1216c7afd0dbSLois Curfman McInnes + snes - the SNES context 1217c7afd0dbSLois Curfman McInnes . func - monitoring routine 1218b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1219b3006f0bSLois Curfman McInnes monitor routine (use PETSC_NULL if no context is desitre) 1220b3006f0bSLois Curfman McInnes - monitordestroy - [optional] routine that frees monitor context 1221b3006f0bSLois Curfman McInnes (may be PETSC_NULL) 12229b94acceSBarry Smith 1223c7afd0dbSLois Curfman McInnes Calling sequence of func: 1224329f5518SBarry Smith $ int func(SNES snes,int its, PetscReal norm,void *mctx) 1225c7afd0dbSLois Curfman McInnes 1226c7afd0dbSLois Curfman McInnes + snes - the SNES context 1227c7afd0dbSLois Curfman McInnes . its - iteration number 1228c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 122940a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 12309b94acceSBarry Smith 12319665c990SLois Curfman McInnes Options Database Keys: 1232c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1233c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1234c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1235c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1236c7afd0dbSLois Curfman McInnes been hardwired into a code by 1237c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1238c7afd0dbSLois Curfman McInnes does not cancel those set via 1239c7afd0dbSLois Curfman McInnes the options database. 12409665c990SLois Curfman McInnes 1241639f9d9dSBarry Smith Notes: 12426bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 12436bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 12446bc08f3fSLois Curfman McInnes order in which they were set. 1245639f9d9dSBarry Smith 124636851e7fSLois Curfman McInnes Level: intermediate 124736851e7fSLois Curfman McInnes 12489b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 12499b94acceSBarry Smith 12505cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 12519b94acceSBarry Smith @*/ 1252329f5518SBarry Smith int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *)) 12539b94acceSBarry Smith { 12543a40ed3dSBarry Smith PetscFunctionBegin; 1255184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1256639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 125729bbc08cSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1258639f9d9dSBarry Smith } 1259639f9d9dSBarry Smith 1260639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1261b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1262639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 12633a40ed3dSBarry Smith PetscFunctionReturn(0); 12649b94acceSBarry Smith } 12659b94acceSBarry Smith 12664a2ae208SSatish Balay #undef __FUNCT__ 12674a2ae208SSatish Balay #define __FUNCT__ "SNESClearMonitor" 12685cd90555SBarry Smith /*@C 12695cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 12705cd90555SBarry Smith 1271c7afd0dbSLois Curfman McInnes Collective on SNES 1272c7afd0dbSLois Curfman McInnes 12735cd90555SBarry Smith Input Parameters: 12745cd90555SBarry Smith . snes - the SNES context 12755cd90555SBarry Smith 12761a480d89SAdministrator Options Database Key: 1277c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1278c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1279c7afd0dbSLois Curfman McInnes set via the options database 12805cd90555SBarry Smith 12815cd90555SBarry Smith Notes: 12825cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 12835cd90555SBarry Smith 128436851e7fSLois Curfman McInnes Level: intermediate 128536851e7fSLois Curfman McInnes 12865cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 12875cd90555SBarry Smith 12885cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 12895cd90555SBarry Smith @*/ 12905cd90555SBarry Smith int SNESClearMonitor(SNES snes) 12915cd90555SBarry Smith { 12925cd90555SBarry Smith PetscFunctionBegin; 1293184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12945cd90555SBarry Smith snes->numbermonitors = 0; 12955cd90555SBarry Smith PetscFunctionReturn(0); 12965cd90555SBarry Smith } 12975cd90555SBarry Smith 12984a2ae208SSatish Balay #undef __FUNCT__ 12994a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceTest" 13009b94acceSBarry Smith /*@C 13019b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13029b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13039b94acceSBarry Smith 1304fee21e36SBarry Smith Collective on SNES 1305fee21e36SBarry Smith 1306c7afd0dbSLois Curfman McInnes Input Parameters: 1307c7afd0dbSLois Curfman McInnes + snes - the SNES context 1308c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1309c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1310c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 13119b94acceSBarry Smith 1312c7afd0dbSLois Curfman McInnes Calling sequence of func: 1313329f5518SBarry Smith $ int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1314c7afd0dbSLois Curfman McInnes 1315c7afd0dbSLois Curfman McInnes + snes - the SNES context 1316c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1317184914b5SBarry Smith . reason - reason for convergence/divergence 1318c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 13194b27c08aSLois Curfman McInnes . gnorm - 2-norm of current step 13204b27c08aSLois Curfman McInnes - f - 2-norm of function 13219b94acceSBarry Smith 132236851e7fSLois Curfman McInnes Level: advanced 132336851e7fSLois Curfman McInnes 13249b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13259b94acceSBarry Smith 13264b27c08aSLois Curfman McInnes .seealso: SNESConverged_LS(), SNESConverged_TR() 13279b94acceSBarry Smith @*/ 1328329f5518SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 13299b94acceSBarry Smith { 13303a40ed3dSBarry Smith PetscFunctionBegin; 1331184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13329b94acceSBarry Smith (snes)->converged = func; 13339b94acceSBarry Smith (snes)->cnvP = cctx; 13343a40ed3dSBarry Smith PetscFunctionReturn(0); 13359b94acceSBarry Smith } 13369b94acceSBarry Smith 13374a2ae208SSatish Balay #undef __FUNCT__ 13384a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergedReason" 1339184914b5SBarry Smith /*@C 1340184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1341184914b5SBarry Smith 1342184914b5SBarry Smith Not Collective 1343184914b5SBarry Smith 1344184914b5SBarry Smith Input Parameter: 1345184914b5SBarry Smith . snes - the SNES context 1346184914b5SBarry Smith 1347184914b5SBarry Smith Output Parameter: 1348e090d566SSatish Balay . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1349184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1350184914b5SBarry Smith 1351184914b5SBarry Smith Level: intermediate 1352184914b5SBarry Smith 1353184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1354184914b5SBarry Smith 1355184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1356184914b5SBarry Smith 13574b27c08aSLois Curfman McInnes .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1358184914b5SBarry Smith @*/ 1359184914b5SBarry Smith int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1360184914b5SBarry Smith { 1361184914b5SBarry Smith PetscFunctionBegin; 1362184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1363184914b5SBarry Smith *reason = snes->reason; 1364184914b5SBarry Smith PetscFunctionReturn(0); 1365184914b5SBarry Smith } 1366184914b5SBarry Smith 13674a2ae208SSatish Balay #undef __FUNCT__ 13684a2ae208SSatish Balay #define __FUNCT__ "SNESSetConvergenceHistory" 1369c9005455SLois Curfman McInnes /*@ 1370c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1371c9005455SLois Curfman McInnes 1372fee21e36SBarry Smith Collective on SNES 1373fee21e36SBarry Smith 1374c7afd0dbSLois Curfman McInnes Input Parameters: 1375c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1376c7afd0dbSLois Curfman McInnes . a - array to hold history 1377cd5578b5SBarry Smith . its - integer array holds the number of linear iterations for each solve. 1378758f92a0SBarry Smith . na - size of a and its 137964731454SLois Curfman McInnes - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1380758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1381c7afd0dbSLois Curfman McInnes 1382c9005455SLois Curfman McInnes Notes: 13834b27c08aSLois Curfman McInnes If set, this array will contain the function norms computed 1384c9005455SLois Curfman McInnes at each step. 1385c9005455SLois Curfman McInnes 1386c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1387c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1388c9005455SLois Curfman McInnes during the section of code that is being timed. 1389c9005455SLois Curfman McInnes 139036851e7fSLois Curfman McInnes Level: intermediate 139136851e7fSLois Curfman McInnes 1392c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1393758f92a0SBarry Smith 139408405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1395758f92a0SBarry Smith 1396c9005455SLois Curfman McInnes @*/ 1397329f5518SBarry Smith int SNESSetConvergenceHistory(SNES snes,PetscReal *a,int *its,int na,PetscTruth reset) 1398c9005455SLois Curfman McInnes { 13993a40ed3dSBarry Smith PetscFunctionBegin; 1400c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1401c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1402c9005455SLois Curfman McInnes snes->conv_hist = a; 1403758f92a0SBarry Smith snes->conv_hist_its = its; 1404758f92a0SBarry Smith snes->conv_hist_max = na; 1405758f92a0SBarry Smith snes->conv_hist_reset = reset; 1406758f92a0SBarry Smith PetscFunctionReturn(0); 1407758f92a0SBarry Smith } 1408758f92a0SBarry Smith 14094a2ae208SSatish Balay #undef __FUNCT__ 14104a2ae208SSatish Balay #define __FUNCT__ "SNESGetConvergenceHistory" 14110c4c9dddSBarry Smith /*@C 1412758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1413758f92a0SBarry Smith 1414758f92a0SBarry Smith Collective on SNES 1415758f92a0SBarry Smith 1416758f92a0SBarry Smith Input Parameter: 1417758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1418758f92a0SBarry Smith 1419758f92a0SBarry Smith Output Parameters: 1420758f92a0SBarry Smith . a - array to hold history 1421758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1422758f92a0SBarry Smith negative if not converged) for each solve. 1423758f92a0SBarry Smith - na - size of a and its 1424758f92a0SBarry Smith 1425758f92a0SBarry Smith Notes: 1426758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1427758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1428758f92a0SBarry Smith 1429758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1430758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1431758f92a0SBarry Smith during the section of code that is being timed. 1432758f92a0SBarry Smith 1433758f92a0SBarry Smith Level: intermediate 1434758f92a0SBarry Smith 1435758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1436758f92a0SBarry Smith 1437758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1438758f92a0SBarry Smith 1439758f92a0SBarry Smith @*/ 1440329f5518SBarry Smith int SNESGetConvergenceHistory(SNES snes,PetscReal **a,int **its,int *na) 1441758f92a0SBarry Smith { 1442758f92a0SBarry Smith PetscFunctionBegin; 1443758f92a0SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1444758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1445758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1446758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 14473a40ed3dSBarry Smith PetscFunctionReturn(0); 1448c9005455SLois Curfman McInnes } 1449c9005455SLois Curfman McInnes 1450e74ef692SMatthew Knepley #undef __FUNCT__ 1451e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetRhsBC" 145276b2cf59SMatthew Knepley /*@ 145376b2cf59SMatthew Knepley SNESSetRhsBC - Sets the function which applies boundary conditions 145476b2cf59SMatthew Knepley to the Rhs of each system. 145576b2cf59SMatthew Knepley 145676b2cf59SMatthew Knepley Collective on SNES 145776b2cf59SMatthew Knepley 145876b2cf59SMatthew Knepley Input Parameters: 145976b2cf59SMatthew Knepley . snes - The nonlinear solver context 146076b2cf59SMatthew Knepley . func - The function 146176b2cf59SMatthew Knepley 146276b2cf59SMatthew Knepley Calling sequence of func: 146376b2cf59SMatthew Knepley . func (SNES snes, Vec rhs, void *ctx); 146476b2cf59SMatthew Knepley 146576b2cf59SMatthew Knepley . rhs - The current rhs vector 146676b2cf59SMatthew Knepley . ctx - The user-context 146776b2cf59SMatthew Knepley 146876b2cf59SMatthew Knepley Level: intermediate 146976b2cf59SMatthew Knepley 147076b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 147176b2cf59SMatthew Knepley .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate() 147276b2cf59SMatthew Knepley @*/ 147376b2cf59SMatthew Knepley int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *)) 147476b2cf59SMatthew Knepley { 147576b2cf59SMatthew Knepley PetscFunctionBegin; 147676b2cf59SMatthew Knepley PetscValidHeaderSpecific(snes, SNES_COOKIE); 147776b2cf59SMatthew Knepley snes->applyrhsbc = func; 147876b2cf59SMatthew Knepley PetscFunctionReturn(0); 147976b2cf59SMatthew Knepley } 148076b2cf59SMatthew Knepley 1481e74ef692SMatthew Knepley #undef __FUNCT__ 1482e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultRhsBC" 148376b2cf59SMatthew Knepley /*@ 148476b2cf59SMatthew Knepley SNESDefaultRhsBC - The default boundary condition function which does nothing. 148576b2cf59SMatthew Knepley 148676b2cf59SMatthew Knepley Not collective 148776b2cf59SMatthew Knepley 148876b2cf59SMatthew Knepley Input Parameters: 148976b2cf59SMatthew Knepley . snes - The nonlinear solver context 149076b2cf59SMatthew Knepley . rhs - The Rhs 149176b2cf59SMatthew Knepley . ctx - The user-context 149276b2cf59SMatthew Knepley 149376b2cf59SMatthew Knepley Level: beginner 149476b2cf59SMatthew Knepley 149576b2cf59SMatthew Knepley .keywords: SNES, Rhs, boundary conditions 149676b2cf59SMatthew Knepley .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate() 149776b2cf59SMatthew Knepley @*/ 149876b2cf59SMatthew Knepley int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx) 149976b2cf59SMatthew Knepley { 150076b2cf59SMatthew Knepley PetscFunctionBegin; 150176b2cf59SMatthew Knepley PetscFunctionReturn(0); 150276b2cf59SMatthew Knepley } 150376b2cf59SMatthew Knepley 1504e74ef692SMatthew Knepley #undef __FUNCT__ 1505e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetSolutionBC" 150676b2cf59SMatthew Knepley /*@ 150776b2cf59SMatthew Knepley SNESSetSolutionBC - Sets the function which applies boundary conditions 150876b2cf59SMatthew Knepley to the solution of each system. 150976b2cf59SMatthew Knepley 151076b2cf59SMatthew Knepley Collective on SNES 151176b2cf59SMatthew Knepley 151276b2cf59SMatthew Knepley Input Parameters: 151376b2cf59SMatthew Knepley . snes - The nonlinear solver context 151476b2cf59SMatthew Knepley . func - The function 151576b2cf59SMatthew Knepley 151676b2cf59SMatthew Knepley Calling sequence of func: 15179d06fb6bSMatthew Knepley . func (SNES snes, Vec rsol, void *ctx); 151876b2cf59SMatthew Knepley 151976b2cf59SMatthew Knepley . sol - The current solution vector 152076b2cf59SMatthew Knepley . ctx - The user-context 152176b2cf59SMatthew Knepley 152276b2cf59SMatthew Knepley Level: intermediate 152376b2cf59SMatthew Knepley 152476b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 152576b2cf59SMatthew Knepley .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate() 152676b2cf59SMatthew Knepley @*/ 152776b2cf59SMatthew Knepley int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *)) 152876b2cf59SMatthew Knepley { 152976b2cf59SMatthew Knepley PetscFunctionBegin; 153076b2cf59SMatthew Knepley PetscValidHeaderSpecific(snes, SNES_COOKIE); 153176b2cf59SMatthew Knepley snes->applysolbc = func; 153276b2cf59SMatthew Knepley PetscFunctionReturn(0); 153376b2cf59SMatthew Knepley } 153476b2cf59SMatthew Knepley 1535e74ef692SMatthew Knepley #undef __FUNCT__ 1536e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultSolutionBC" 153776b2cf59SMatthew Knepley /*@ 153876b2cf59SMatthew Knepley SNESDefaultSolutionBC - The default boundary condition function which does nothing. 153976b2cf59SMatthew Knepley 154076b2cf59SMatthew Knepley Not collective 154176b2cf59SMatthew Knepley 154276b2cf59SMatthew Knepley Input Parameters: 154376b2cf59SMatthew Knepley . snes - The nonlinear solver context 154476b2cf59SMatthew Knepley . sol - The solution 154576b2cf59SMatthew Knepley . ctx - The user-context 154676b2cf59SMatthew Knepley 154776b2cf59SMatthew Knepley Level: beginner 154876b2cf59SMatthew Knepley 154976b2cf59SMatthew Knepley .keywords: SNES, solution, boundary conditions 155076b2cf59SMatthew Knepley .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate() 155176b2cf59SMatthew Knepley @*/ 155276b2cf59SMatthew Knepley int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx) 155376b2cf59SMatthew Knepley { 155476b2cf59SMatthew Knepley PetscFunctionBegin; 155576b2cf59SMatthew Knepley PetscFunctionReturn(0); 155676b2cf59SMatthew Knepley } 155776b2cf59SMatthew Knepley 1558e74ef692SMatthew Knepley #undef __FUNCT__ 1559e74ef692SMatthew Knepley #define __FUNCT__ "SNESSetUpdate" 156076b2cf59SMatthew Knepley /*@ 156176b2cf59SMatthew Knepley SNESSetUpdate - Sets the general-purpose update function called 156276b2cf59SMatthew Knepley at the beginning of every step of the iteration. 156376b2cf59SMatthew Knepley 156476b2cf59SMatthew Knepley Collective on SNES 156576b2cf59SMatthew Knepley 156676b2cf59SMatthew Knepley Input Parameters: 156776b2cf59SMatthew Knepley . snes - The nonlinear solver context 156876b2cf59SMatthew Knepley . func - The function 156976b2cf59SMatthew Knepley 157076b2cf59SMatthew Knepley Calling sequence of func: 157176b2cf59SMatthew Knepley . func (TS ts, int step); 157276b2cf59SMatthew Knepley 157376b2cf59SMatthew Knepley . step - The current step of the iteration 157476b2cf59SMatthew Knepley 157576b2cf59SMatthew Knepley Level: intermediate 157676b2cf59SMatthew Knepley 157776b2cf59SMatthew Knepley .keywords: SNES, update 157876b2cf59SMatthew Knepley .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 157976b2cf59SMatthew Knepley @*/ 158076b2cf59SMatthew Knepley int SNESSetUpdate(SNES snes, int (*func)(SNES, int)) 158176b2cf59SMatthew Knepley { 158276b2cf59SMatthew Knepley PetscFunctionBegin; 158376b2cf59SMatthew Knepley PetscValidHeaderSpecific(snes, SNES_COOKIE); 158476b2cf59SMatthew Knepley snes->update = func; 158576b2cf59SMatthew Knepley PetscFunctionReturn(0); 158676b2cf59SMatthew Knepley } 158776b2cf59SMatthew Knepley 1588e74ef692SMatthew Knepley #undef __FUNCT__ 1589e74ef692SMatthew Knepley #define __FUNCT__ "SNESDefaultUpdate" 159076b2cf59SMatthew Knepley /*@ 159176b2cf59SMatthew Knepley SNESDefaultUpdate - The default update function which does nothing. 159276b2cf59SMatthew Knepley 159376b2cf59SMatthew Knepley Not collective 159476b2cf59SMatthew Knepley 159576b2cf59SMatthew Knepley Input Parameters: 159676b2cf59SMatthew Knepley . snes - The nonlinear solver context 159776b2cf59SMatthew Knepley . step - The current step of the iteration 159876b2cf59SMatthew Knepley 1599205452f4SMatthew Knepley Level: intermediate 1600205452f4SMatthew Knepley 160176b2cf59SMatthew Knepley .keywords: SNES, update 160276b2cf59SMatthew Knepley .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 160376b2cf59SMatthew Knepley @*/ 160476b2cf59SMatthew Knepley int SNESDefaultUpdate(SNES snes, int step) 160576b2cf59SMatthew Knepley { 160676b2cf59SMatthew Knepley PetscFunctionBegin; 160776b2cf59SMatthew Knepley PetscFunctionReturn(0); 160876b2cf59SMatthew Knepley } 160976b2cf59SMatthew Knepley 16104a2ae208SSatish Balay #undef __FUNCT__ 16114a2ae208SSatish Balay #define __FUNCT__ "SNESScaleStep_Private" 16129b94acceSBarry Smith /* 16139b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 16149b94acceSBarry Smith positive parameter delta. 16159b94acceSBarry Smith 16169b94acceSBarry Smith Input Parameters: 1617c7afd0dbSLois Curfman McInnes + snes - the SNES context 16189b94acceSBarry Smith . y - approximate solution of linear system 16199b94acceSBarry Smith . fnorm - 2-norm of current function 1620c7afd0dbSLois Curfman McInnes - delta - trust region size 16219b94acceSBarry Smith 16229b94acceSBarry Smith Output Parameters: 1623c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 16249b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 16259b94acceSBarry Smith region, and exceeds zero otherwise. 1626c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 16279b94acceSBarry Smith 16289b94acceSBarry Smith Note: 16294b27c08aSLois Curfman McInnes For non-trust region methods such as SNESLS, the parameter delta 16309b94acceSBarry Smith is set to be the maximum allowable step size. 16319b94acceSBarry Smith 16329b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16339b94acceSBarry Smith */ 1634064f8208SBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 16359b94acceSBarry Smith { 1636064f8208SBarry Smith PetscReal nrm; 1637ea709b57SSatish Balay PetscScalar cnorm; 16383a40ed3dSBarry Smith int ierr; 16393a40ed3dSBarry Smith 16403a40ed3dSBarry Smith PetscFunctionBegin; 1641184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1642184914b5SBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE); 1643184914b5SBarry Smith PetscCheckSameComm(snes,y); 1644184914b5SBarry Smith 1645064f8208SBarry Smith ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1646064f8208SBarry Smith if (nrm > *delta) { 1647064f8208SBarry Smith nrm = *delta/nrm; 1648064f8208SBarry Smith *gpnorm = (1.0 - nrm)*(*fnorm); 1649064f8208SBarry Smith cnorm = nrm; 1650329f5518SBarry Smith ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 16519b94acceSBarry Smith *ynorm = *delta; 16529b94acceSBarry Smith } else { 16539b94acceSBarry Smith *gpnorm = 0.0; 1654064f8208SBarry Smith *ynorm = nrm; 16559b94acceSBarry Smith } 16563a40ed3dSBarry Smith PetscFunctionReturn(0); 16579b94acceSBarry Smith } 16589b94acceSBarry Smith 16594a2ae208SSatish Balay #undef __FUNCT__ 16604a2ae208SSatish Balay #define __FUNCT__ "SNESSolve" 16619b94acceSBarry Smith /*@ 16629b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 166363a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 16649b94acceSBarry Smith 1665c7afd0dbSLois Curfman McInnes Collective on SNES 1666c7afd0dbSLois Curfman McInnes 1667b2002411SLois Curfman McInnes Input Parameters: 1668c7afd0dbSLois Curfman McInnes + snes - the SNES context 1669c7afd0dbSLois Curfman McInnes - x - the solution vector 16709b94acceSBarry Smith 16719b94acceSBarry Smith Output Parameter: 1672b2002411SLois Curfman McInnes . its - number of iterations until termination 16739b94acceSBarry Smith 1674b2002411SLois Curfman McInnes Notes: 16758ddd3da0SLois Curfman McInnes The user should initialize the vector,x, with the initial guess 16768ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 16778ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 16788ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 16798ddd3da0SLois Curfman McInnes 168036851e7fSLois Curfman McInnes Level: beginner 168136851e7fSLois Curfman McInnes 16829b94acceSBarry Smith .keywords: SNES, nonlinear, solve 16839b94acceSBarry Smith 168463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 16859b94acceSBarry Smith @*/ 16868ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 16879b94acceSBarry Smith { 1688f1af5d2fSBarry Smith int ierr; 1689f1af5d2fSBarry Smith PetscTruth flg; 1690052efed2SBarry Smith 16913a40ed3dSBarry Smith PetscFunctionBegin; 169277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1693184914b5SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 1694184914b5SBarry Smith PetscCheckSameComm(snes,x); 169574679c65SBarry Smith PetscValidIntPointer(its); 169629bbc08cSBarry Smith if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 169774637425SBarry Smith 169882bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1699c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 1700758f92a0SBarry Smith if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 1701d5ba7fb7SMatthew Knepley ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 170250ffb88aSMatthew Knepley snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 17039b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr); 1704d5ba7fb7SMatthew Knepley ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1705b0a32e0cSBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 17068b179ff8SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1707da9b6338SBarry Smith ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1708da9b6338SBarry Smith if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 17093a40ed3dSBarry Smith PetscFunctionReturn(0); 17109b94acceSBarry Smith } 17119b94acceSBarry Smith 17129b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 17139b94acceSBarry Smith 17144a2ae208SSatish Balay #undef __FUNCT__ 17154a2ae208SSatish Balay #define __FUNCT__ "SNESSetType" 171682bf6240SBarry Smith /*@C 17174b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 17189b94acceSBarry Smith 1719fee21e36SBarry Smith Collective on SNES 1720fee21e36SBarry Smith 1721c7afd0dbSLois Curfman McInnes Input Parameters: 1722c7afd0dbSLois Curfman McInnes + snes - the SNES context 1723454a90a3SBarry Smith - type - a known method 1724c7afd0dbSLois Curfman McInnes 1725c7afd0dbSLois Curfman McInnes Options Database Key: 1726454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 1727c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1728ae12b187SLois Curfman McInnes 17299b94acceSBarry Smith Notes: 1730e090d566SSatish Balay See "petsc/include/petscsnes.h" for available methods (for instance) 17314b27c08aSLois Curfman McInnes + SNESLS - Newton's method with line search 1732c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17334b27c08aSLois Curfman McInnes . SNESTR - Newton's method with trust region 1734c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 17359b94acceSBarry Smith 1736ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1737ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1738ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1739ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1740ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1741ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1742ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1743ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1744ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1745b0a32e0cSBarry Smith appropriate method. 174636851e7fSLois Curfman McInnes 174736851e7fSLois Curfman McInnes Level: intermediate 1748a703fe33SLois Curfman McInnes 1749454a90a3SBarry Smith .keywords: SNES, set, type 1750435da068SBarry Smith 1751435da068SBarry Smith .seealso: SNESType, SNESCreate() 1752435da068SBarry Smith 17539b94acceSBarry Smith @*/ 1754454a90a3SBarry Smith int SNESSetType(SNES snes,SNESType type) 17559b94acceSBarry Smith { 17560f5bd95cSBarry Smith int ierr,(*r)(SNES); 17576831982aSBarry Smith PetscTruth match; 17583a40ed3dSBarry Smith 17593a40ed3dSBarry Smith PetscFunctionBegin; 176077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17610f5bd95cSBarry Smith PetscValidCharPointer(type); 176282bf6240SBarry Smith 17636831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 17640f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 176582bf6240SBarry Smith 176682bf6240SBarry Smith if (snes->setupcalled) { 1767e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 176882bf6240SBarry Smith snes->data = 0; 176982bf6240SBarry Smith } 17707d1a2b2bSBarry Smith 17719b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 177282bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 177382bf6240SBarry Smith 1774b9617806SBarry Smith ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 177582bf6240SBarry Smith 177629bbc08cSBarry Smith if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type); 17771548b14aSBarry Smith 1778606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 177982bf6240SBarry Smith snes->data = 0; 17803a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 1781454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 178282bf6240SBarry Smith 17833a40ed3dSBarry Smith PetscFunctionReturn(0); 17849b94acceSBarry Smith } 17859b94acceSBarry Smith 1786a847f771SSatish Balay 17879b94acceSBarry Smith /* --------------------------------------------------------------------- */ 17884a2ae208SSatish Balay #undef __FUNCT__ 17894a2ae208SSatish Balay #define __FUNCT__ "SNESRegisterDestroy" 17909b94acceSBarry Smith /*@C 17919b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1792f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 17939b94acceSBarry Smith 1794fee21e36SBarry Smith Not Collective 1795fee21e36SBarry Smith 179636851e7fSLois Curfman McInnes Level: advanced 179736851e7fSLois Curfman McInnes 17989b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 17999b94acceSBarry Smith 18009b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 18019b94acceSBarry Smith @*/ 1802cf256101SBarry Smith int SNESRegisterDestroy(void) 18039b94acceSBarry Smith { 180482bf6240SBarry Smith int ierr; 180582bf6240SBarry Smith 18063a40ed3dSBarry Smith PetscFunctionBegin; 180782bf6240SBarry Smith if (SNESList) { 1808b0a32e0cSBarry Smith ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 180982bf6240SBarry Smith SNESList = 0; 18109b94acceSBarry Smith } 18114c49b128SBarry Smith SNESRegisterAllCalled = PETSC_FALSE; 18123a40ed3dSBarry Smith PetscFunctionReturn(0); 18139b94acceSBarry Smith } 18149b94acceSBarry Smith 18154a2ae208SSatish Balay #undef __FUNCT__ 18164a2ae208SSatish Balay #define __FUNCT__ "SNESGetType" 18179b94acceSBarry Smith /*@C 18189a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 18199b94acceSBarry Smith 1820c7afd0dbSLois Curfman McInnes Not Collective 1821c7afd0dbSLois Curfman McInnes 18229b94acceSBarry Smith Input Parameter: 18234b0e389bSBarry Smith . snes - nonlinear solver context 18249b94acceSBarry Smith 18259b94acceSBarry Smith Output Parameter: 18263a7fca6bSBarry Smith . type - SNES method (a character string) 18279b94acceSBarry Smith 182836851e7fSLois Curfman McInnes Level: intermediate 182936851e7fSLois Curfman McInnes 1830454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 18319b94acceSBarry Smith @*/ 1832454a90a3SBarry Smith int SNESGetType(SNES snes,SNESType *type) 18339b94acceSBarry Smith { 18343a40ed3dSBarry Smith PetscFunctionBegin; 1835184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1836454a90a3SBarry Smith *type = snes->type_name; 18373a40ed3dSBarry Smith PetscFunctionReturn(0); 18389b94acceSBarry Smith } 18399b94acceSBarry Smith 18404a2ae208SSatish Balay #undef __FUNCT__ 18414a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolution" 18429b94acceSBarry Smith /*@C 18439b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 18449b94acceSBarry Smith stored. 18459b94acceSBarry Smith 1846c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1847c7afd0dbSLois Curfman McInnes 18489b94acceSBarry Smith Input Parameter: 18499b94acceSBarry Smith . snes - the SNES context 18509b94acceSBarry Smith 18519b94acceSBarry Smith Output Parameter: 18529b94acceSBarry Smith . x - the solution 18539b94acceSBarry Smith 185436851e7fSLois Curfman McInnes Level: advanced 185536851e7fSLois Curfman McInnes 18569b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 18579b94acceSBarry Smith 18584b27c08aSLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetSolutionUpdate() 18599b94acceSBarry Smith @*/ 18609b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 18619b94acceSBarry Smith { 18623a40ed3dSBarry Smith PetscFunctionBegin; 186377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18649b94acceSBarry Smith *x = snes->vec_sol_always; 18653a40ed3dSBarry Smith PetscFunctionReturn(0); 18669b94acceSBarry Smith } 18679b94acceSBarry Smith 18684a2ae208SSatish Balay #undef __FUNCT__ 18694a2ae208SSatish Balay #define __FUNCT__ "SNESGetSolutionUpdate" 18709b94acceSBarry Smith /*@C 18719b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 18729b94acceSBarry Smith stored. 18739b94acceSBarry Smith 1874c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1875c7afd0dbSLois Curfman McInnes 18769b94acceSBarry Smith Input Parameter: 18779b94acceSBarry Smith . snes - the SNES context 18789b94acceSBarry Smith 18799b94acceSBarry Smith Output Parameter: 18809b94acceSBarry Smith . x - the solution update 18819b94acceSBarry Smith 188236851e7fSLois Curfman McInnes Level: advanced 188336851e7fSLois Curfman McInnes 18849b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 18859b94acceSBarry Smith 18869b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 18879b94acceSBarry Smith @*/ 18889b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 18899b94acceSBarry Smith { 18903a40ed3dSBarry Smith PetscFunctionBegin; 189177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18929b94acceSBarry Smith *x = snes->vec_sol_update_always; 18933a40ed3dSBarry Smith PetscFunctionReturn(0); 18949b94acceSBarry Smith } 18959b94acceSBarry Smith 18964a2ae208SSatish Balay #undef __FUNCT__ 18974a2ae208SSatish Balay #define __FUNCT__ "SNESGetFunction" 18989b94acceSBarry Smith /*@C 18993638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 19009b94acceSBarry Smith 1901c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1902c7afd0dbSLois Curfman McInnes 19039b94acceSBarry Smith Input Parameter: 19049b94acceSBarry Smith . snes - the SNES context 19059b94acceSBarry Smith 19069b94acceSBarry Smith Output Parameter: 19077bf4e008SBarry Smith + r - the function (or PETSC_NULL) 190800036973SBarry Smith . ctx - the function context (or PETSC_NULL) 190900036973SBarry Smith - func - the function (or PETSC_NULL) 19109b94acceSBarry Smith 191136851e7fSLois Curfman McInnes Level: advanced 191236851e7fSLois Curfman McInnes 1913a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 19149b94acceSBarry Smith 19154b27c08aSLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution() 19169b94acceSBarry Smith @*/ 191700036973SBarry Smith int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*)) 19189b94acceSBarry Smith { 19193a40ed3dSBarry Smith PetscFunctionBegin; 192077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 19217bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 19227bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 192300036973SBarry Smith if (func) *func = snes->computefunction; 19243a40ed3dSBarry Smith PetscFunctionReturn(0); 19259b94acceSBarry Smith } 19269b94acceSBarry Smith 19274a2ae208SSatish Balay #undef __FUNCT__ 19284a2ae208SSatish Balay #define __FUNCT__ "SNESSetOptionsPrefix" 19293c7409f5SSatish Balay /*@C 19303c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1931d850072dSLois Curfman McInnes SNES options in the database. 19323c7409f5SSatish Balay 1933fee21e36SBarry Smith Collective on SNES 1934fee21e36SBarry Smith 1935c7afd0dbSLois Curfman McInnes Input Parameter: 1936c7afd0dbSLois Curfman McInnes + snes - the SNES context 1937c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1938c7afd0dbSLois Curfman McInnes 1939d850072dSLois Curfman McInnes Notes: 1940a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1941c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1942d850072dSLois Curfman McInnes 194336851e7fSLois Curfman McInnes Level: advanced 194436851e7fSLois Curfman McInnes 19453c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1946a86d99e1SLois Curfman McInnes 1947a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19483c7409f5SSatish Balay @*/ 19493c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 19503c7409f5SSatish Balay { 19513c7409f5SSatish Balay int ierr; 19523c7409f5SSatish Balay 19533a40ed3dSBarry Smith PetscFunctionBegin; 195477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1955639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 19563c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19573a40ed3dSBarry Smith PetscFunctionReturn(0); 19583c7409f5SSatish Balay } 19593c7409f5SSatish Balay 19604a2ae208SSatish Balay #undef __FUNCT__ 19614a2ae208SSatish Balay #define __FUNCT__ "SNESAppendOptionsPrefix" 19623c7409f5SSatish Balay /*@C 1963f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1964d850072dSLois Curfman McInnes SNES options in the database. 19653c7409f5SSatish Balay 1966fee21e36SBarry Smith Collective on SNES 1967fee21e36SBarry Smith 1968c7afd0dbSLois Curfman McInnes Input Parameters: 1969c7afd0dbSLois Curfman McInnes + snes - the SNES context 1970c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1971c7afd0dbSLois Curfman McInnes 1972d850072dSLois Curfman McInnes Notes: 1973a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1974c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1975d850072dSLois Curfman McInnes 197636851e7fSLois Curfman McInnes Level: advanced 197736851e7fSLois Curfman McInnes 19783c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1979a86d99e1SLois Curfman McInnes 1980a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 19813c7409f5SSatish Balay @*/ 19823c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 19833c7409f5SSatish Balay { 19843c7409f5SSatish Balay int ierr; 19853c7409f5SSatish Balay 19863a40ed3dSBarry Smith PetscFunctionBegin; 198777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1988639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 19893c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19903a40ed3dSBarry Smith PetscFunctionReturn(0); 19913c7409f5SSatish Balay } 19923c7409f5SSatish Balay 19934a2ae208SSatish Balay #undef __FUNCT__ 19944a2ae208SSatish Balay #define __FUNCT__ "SNESGetOptionsPrefix" 19959ab63eb5SSatish Balay /*@C 19963c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 19973c7409f5SSatish Balay SNES options in the database. 19983c7409f5SSatish Balay 1999c7afd0dbSLois Curfman McInnes Not Collective 2000c7afd0dbSLois Curfman McInnes 20013c7409f5SSatish Balay Input Parameter: 20023c7409f5SSatish Balay . snes - the SNES context 20033c7409f5SSatish Balay 20043c7409f5SSatish Balay Output Parameter: 20053c7409f5SSatish Balay . prefix - pointer to the prefix string used 20063c7409f5SSatish Balay 20079ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 20089ab63eb5SSatish Balay sufficient length to hold the prefix. 20099ab63eb5SSatish Balay 201036851e7fSLois Curfman McInnes Level: advanced 201136851e7fSLois Curfman McInnes 20123c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2013a86d99e1SLois Curfman McInnes 2014a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20153c7409f5SSatish Balay @*/ 20163c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 20173c7409f5SSatish Balay { 20183c7409f5SSatish Balay int ierr; 20193c7409f5SSatish Balay 20203a40ed3dSBarry Smith PetscFunctionBegin; 202177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2022639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 20233a40ed3dSBarry Smith PetscFunctionReturn(0); 20243c7409f5SSatish Balay } 20253c7409f5SSatish Balay 2026acb85ae6SSatish Balay /*MC 2027f1af5d2fSBarry Smith SNESRegisterDynamic - Adds a method to the nonlinear solver package. 20289b94acceSBarry Smith 2029b2002411SLois Curfman McInnes Synopsis: 20303a7fca6bSBarry Smith int SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES)) 20319b94acceSBarry Smith 20328d76a1e5SLois Curfman McInnes Not collective 20338d76a1e5SLois Curfman McInnes 2034b2002411SLois Curfman McInnes Input Parameters: 2035c7afd0dbSLois Curfman McInnes + name_solver - name of a new user-defined solver 2036b2002411SLois Curfman McInnes . path - path (either absolute or relative) the library containing this solver 2037b2002411SLois Curfman McInnes . name_create - name of routine to create method context 2038c7afd0dbSLois Curfman McInnes - routine_create - routine to create method context 20399b94acceSBarry Smith 2040b2002411SLois Curfman McInnes Notes: 2041f1af5d2fSBarry Smith SNESRegisterDynamic() may be called multiple times to add several user-defined solvers. 20423a40ed3dSBarry Smith 2043b2002411SLois Curfman McInnes If dynamic libraries are used, then the fourth input argument (routine_create) 2044b2002411SLois Curfman McInnes is ignored. 2045b2002411SLois Curfman McInnes 2046b9617806SBarry Smith Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR}, ${BOPT}, 20475ba88a07SLois Curfman McInnes and others of the form ${any_environmental_variable} occuring in pathname will be 20485ba88a07SLois Curfman McInnes replaced with appropriate values. 20495ba88a07SLois Curfman McInnes 2050b2002411SLois Curfman McInnes Sample usage: 205169225d78SLois Curfman McInnes .vb 2052f1af5d2fSBarry Smith SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a, 2053b2002411SLois Curfman McInnes "MySolverCreate",MySolverCreate); 205469225d78SLois Curfman McInnes .ve 2055b2002411SLois Curfman McInnes 2056b2002411SLois Curfman McInnes Then, your solver can be chosen with the procedural interface via 2057b2002411SLois Curfman McInnes $ SNESSetType(snes,"my_solver") 2058b2002411SLois Curfman McInnes or at runtime via the option 2059b2002411SLois Curfman McInnes $ -snes_type my_solver 2060b2002411SLois Curfman McInnes 206136851e7fSLois Curfman McInnes Level: advanced 206236851e7fSLois Curfman McInnes 2063*3cea93caSBarry Smith Note: If your function is not being put into a shared library then use SNESRegister() instead 2064*3cea93caSBarry Smith 2065b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register 2066b2002411SLois Curfman McInnes 2067b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy() 2068b2002411SLois Curfman McInnes M*/ 2069b2002411SLois Curfman McInnes 20704a2ae208SSatish Balay #undef __FUNCT__ 20714a2ae208SSatish Balay #define __FUNCT__ "SNESRegister" 2072*3cea93caSBarry Smith /*@C 2073*3cea93caSBarry Smith SNESRegister - See SNESRegisterDynamic() 2074*3cea93caSBarry Smith 2075*3cea93caSBarry Smith @*/ 2076f1af5d2fSBarry Smith int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES)) 2077b2002411SLois Curfman McInnes { 2078b2002411SLois Curfman McInnes char fullname[256]; 2079b2002411SLois Curfman McInnes int ierr; 2080b2002411SLois Curfman McInnes 2081b2002411SLois Curfman McInnes PetscFunctionBegin; 2082b0a32e0cSBarry Smith ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2083c134de8dSSatish Balay ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2084b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2085b2002411SLois Curfman McInnes } 2086da9b6338SBarry Smith 2087da9b6338SBarry Smith #undef __FUNCT__ 2088da9b6338SBarry Smith #define __FUNCT__ "SNESTestLocalMin" 2089da9b6338SBarry Smith int SNESTestLocalMin(SNES snes) 2090da9b6338SBarry Smith { 2091da9b6338SBarry Smith int ierr,N,i,j; 2092da9b6338SBarry Smith Vec u,uh,fh; 2093da9b6338SBarry Smith PetscScalar value; 2094da9b6338SBarry Smith PetscReal norm; 2095da9b6338SBarry Smith 2096da9b6338SBarry Smith PetscFunctionBegin; 2097da9b6338SBarry Smith ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2098da9b6338SBarry Smith ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2099da9b6338SBarry Smith ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2100da9b6338SBarry Smith 2101da9b6338SBarry Smith /* currently only works for sequential */ 2102da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2103da9b6338SBarry Smith ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2104da9b6338SBarry Smith for (i=0; i<N; i++) { 2105da9b6338SBarry Smith ierr = VecCopy(u,uh);CHKERRQ(ierr); 2106da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %d\n",i);CHKERRQ(ierr); 2107da9b6338SBarry Smith for (j=-10; j<11; j++) { 2108ccae9161SBarry Smith value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2109da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2110da9b6338SBarry Smith ierr = (*snes->computefunction)(snes,uh,fh,snes->funP);CHKERRQ(ierr); 2111da9b6338SBarry Smith ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2112da9b6338SBarry Smith ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %d %18.16e\n",j,norm);CHKERRQ(ierr); 2113da9b6338SBarry Smith value = -value; 2114da9b6338SBarry Smith ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2115da9b6338SBarry Smith } 2116da9b6338SBarry Smith } 2117da9b6338SBarry Smith ierr = VecDestroy(uh);CHKERRQ(ierr); 2118da9b6338SBarry Smith ierr = VecDestroy(fh);CHKERRQ(ierr); 2119da9b6338SBarry Smith PetscFunctionReturn(0); 2120da9b6338SBarry Smith } 2121