1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*a8c6a408SBarry Smith static char vcid[] = "$Id: snes.c,v 1.134 1997/11/28 16:21:38 bsmith Exp bsmith $"; 39b94acceSBarry Smith #endif 49b94acceSBarry Smith 570f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 6f5eb4b81SSatish Balay #include "src/sys/nreg.h" 79b94acceSBarry Smith #include "pinclude/pviewer.h" 89b94acceSBarry Smith #include <math.h> 99b94acceSBarry Smith 1084cb2905SBarry Smith int SNESRegisterAllCalled = 0; 11ca161407SBarry Smith static NRList *__SNESList = 0; 1284cb2905SBarry Smith 135615d1e5SSatish Balay #undef __FUNC__ 14d4bb536fSBarry Smith #define __FUNC__ "SNESView" 159b94acceSBarry Smith /*@ 169b94acceSBarry Smith SNESView - Prints the SNES data structure. 179b94acceSBarry Smith 189b94acceSBarry Smith Input Parameters: 199b94acceSBarry Smith . SNES - the SNES context 209b94acceSBarry Smith . viewer - visualization context 219b94acceSBarry Smith 229b94acceSBarry Smith Options Database Key: 239b94acceSBarry Smith $ -snes_view : calls SNESView() at end of SNESSolve() 249b94acceSBarry Smith 259b94acceSBarry Smith Notes: 269b94acceSBarry Smith The available visualization contexts include 276d4a8577SBarry Smith $ VIEWER_STDOUT_SELF - standard output (default) 286d4a8577SBarry Smith $ VIEWER_STDOUT_WORLD - synchronized standard 299b94acceSBarry Smith $ output where only the first processor opens 309b94acceSBarry Smith $ the file. All other processors send their 319b94acceSBarry Smith $ data to the first processor to print. 329b94acceSBarry Smith 339b94acceSBarry Smith The user can open alternative vistualization contexts with 34dbb450caSBarry Smith $ ViewerFileOpenASCII() - output to a specified file 359b94acceSBarry Smith 369b94acceSBarry Smith .keywords: SNES, view 379b94acceSBarry Smith 38dbb450caSBarry Smith .seealso: ViewerFileOpenASCII() 399b94acceSBarry Smith @*/ 409b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 419b94acceSBarry Smith { 429b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 439b94acceSBarry Smith FILE *fd; 449b94acceSBarry Smith int ierr; 459b94acceSBarry Smith SLES sles; 469b94acceSBarry Smith char *method; 4719bcc07fSBarry Smith ViewerType vtype; 489b94acceSBarry Smith 493a40ed3dSBarry Smith PetscFunctionBegin; 5074679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5174679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5274679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5374679c65SBarry Smith 5419bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 5519bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 5690ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 5777c4ece6SBarry Smith PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 584b0e389bSBarry Smith SNESGetType(snes,PETSC_NULL,&method); 5977c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," method: %s\n",method); 609b94acceSBarry Smith if (snes->view) (*snes->view)((PetscObject)snes,viewer); 6177c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 629b94acceSBarry Smith " maximum iterations=%d, maximum function evaluations=%d\n", 639b94acceSBarry Smith snes->max_its,snes->max_funcs); 6477c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 659b94acceSBarry Smith " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 669b94acceSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol); 677a00f4a9SLois Curfman McInnes PetscFPrintf(snes->comm,fd, 687a00f4a9SLois Curfman McInnes " total number of linear solver iterations=%d\n",snes->linear_its); 69ae3c334cSLois Curfman McInnes PetscFPrintf(snes->comm,fd, 7068632166SLois Curfman McInnes " total number of function evaluations=%d\n",snes->nfuncs); 719b94acceSBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 7277c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 739b94acceSBarry Smith if (snes->ksp_ewconv) { 749b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 759b94acceSBarry Smith if (kctx) { 7677c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 779b94acceSBarry Smith " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 789b94acceSBarry Smith kctx->version); 7977c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 809b94acceSBarry Smith " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 819b94acceSBarry Smith kctx->rtol_max,kctx->threshold); 8277c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 839b94acceSBarry Smith kctx->gamma,kctx->alpha,kctx->alpha2); 849b94acceSBarry Smith } 859b94acceSBarry Smith } 86c30f285eSLois Curfman McInnes } else if (vtype == STRING_VIEWER) { 87c30f285eSLois Curfman McInnes SNESGetType(snes,PETSC_NULL,&method); 88c30f285eSLois Curfman McInnes ViewerStringSPrintf(viewer," %-3.3s",method); 8919bcc07fSBarry Smith } 909b94acceSBarry Smith SNESGetSLES(snes,&sles); 919b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 923a40ed3dSBarry Smith PetscFunctionReturn(0); 939b94acceSBarry Smith } 949b94acceSBarry Smith 95639f9d9dSBarry Smith /* 96639f9d9dSBarry Smith We retain a list of functions that also take SNES command 97639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 98639f9d9dSBarry Smith */ 99639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 100639f9d9dSBarry Smith static int numberofsetfromoptions; 101639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 102639f9d9dSBarry Smith 1035615d1e5SSatish Balay #undef __FUNC__ 104d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker" 105639f9d9dSBarry Smith /*@ 106639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 107639f9d9dSBarry Smith 108639f9d9dSBarry Smith Input Parameter: 109639f9d9dSBarry Smith . snescheck - function that checks for options 110639f9d9dSBarry Smith 111639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 112639f9d9dSBarry Smith @*/ 113639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 114639f9d9dSBarry Smith { 1153a40ed3dSBarry Smith PetscFunctionBegin; 116639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 117*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed"); 118639f9d9dSBarry Smith } 119639f9d9dSBarry Smith 120639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 1213a40ed3dSBarry Smith PetscFunctionReturn(0); 122639f9d9dSBarry Smith } 123639f9d9dSBarry Smith 1245615d1e5SSatish Balay #undef __FUNC__ 1255615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1269b94acceSBarry Smith /*@ 127682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1289b94acceSBarry Smith 1299b94acceSBarry Smith Input Parameter: 1309b94acceSBarry Smith . snes - the SNES context 1319b94acceSBarry Smith 13211ca99fdSLois Curfman McInnes Notes: 13311ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 13411ca99fdSLois Curfman McInnes the users manual. 13583e2fdc7SBarry Smith 1369b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1379b94acceSBarry Smith 138a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1399b94acceSBarry Smith @*/ 1409b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1419b94acceSBarry Smith { 1424b0e389bSBarry Smith SNESType method; 1439b94acceSBarry Smith double tmp; 1449b94acceSBarry Smith SLES sles; 14581f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1463c7409f5SSatish Balay int version = PETSC_DEFAULT; 1479b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1489b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1499b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1509b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1519b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1529b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1539b94acceSBarry Smith 1543a40ed3dSBarry Smith PetscFunctionBegin; 15581f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 15681f57ec7SBarry Smith 15777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 158*a8c6a408SBarry Smith if (snes->setup_called) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp"); 159ca161407SBarry Smith 160ca161407SBarry Smith if (!__SNESList) {ierr = SNESRegisterAll();CHKERRQ(ierr);} 161ca161407SBarry Smith ierr = NRGetTypeFromOptions(snes->prefix,"-snes_type",__SNESList,&method,&flg);CHKERRQ(ierr); 162052efed2SBarry Smith if (flg) { 163052efed2SBarry Smith ierr = SNESSetType(snes,method); CHKERRQ(ierr); 164ca161407SBarry Smith } else if (!snes->set_method_called) { 1655334005bSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 16640191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 167d64ed03dSBarry Smith } else { 16840191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 1695334005bSBarry Smith } 1705334005bSBarry Smith } 1715334005bSBarry Smith 1723c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 173d64ed03dSBarry Smith if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);} 1743c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 175d64ed03dSBarry Smith if (flg) { 176d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 177d64ed03dSBarry Smith } 1783c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 179d64ed03dSBarry Smith if (flg) { 180d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 181d64ed03dSBarry Smith } 1823c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 183d64ed03dSBarry Smith if (flg) { 184d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 185d64ed03dSBarry Smith } 1863c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 1873c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 188d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 189d64ed03dSBarry Smith if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); } 190d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 191d64ed03dSBarry Smith if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp); CHKERRQ(ierr);} 1923c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 1933c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 194b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 195b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 196b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 197b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 198b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 199b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 200b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 20193c39befSBarry Smith 20293c39befSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr); 20393c39befSBarry Smith if (flg) {snes->converged = 0;} 20493c39befSBarry Smith 2059b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 2069b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 2075bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 2085bbad29bSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);} 2093c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 210639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 2113c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 212d31a3109SLois Curfman McInnes if (flg) { 213639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 214d31a3109SLois Curfman McInnes } 21581f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2163c7409f5SSatish Balay if (flg) { 21717699dbbSLois Curfman McInnes int rank = 0; 218d7e8b826SBarry Smith DrawLG lg; 21917699dbbSLois Curfman McInnes MPI_Initialized(&rank); 22017699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 22117699dbbSLois Curfman McInnes if (!rank) { 22281f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2239b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 224f8c826e1SBarry Smith snes->xmonitor = lg; 2259b94acceSBarry Smith PLogObjectParent(snes,lg); 2269b94acceSBarry Smith } 2279b94acceSBarry Smith } 2283c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2293c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2309b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2319b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 23294a424c1SBarry Smith PLogInfo(snes, 23331615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 23431615d04SLois Curfman McInnes } 23531615d04SLois Curfman McInnes else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 23631615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 23731615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 238d64ed03dSBarry Smith PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2399b94acceSBarry Smith } 240639f9d9dSBarry Smith 241639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 242639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 243639f9d9dSBarry Smith } 244639f9d9dSBarry Smith 2459b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2469b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 2473a40ed3dSBarry Smith if (snes->setfromoptions) { 2483a40ed3dSBarry Smith ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 2493a40ed3dSBarry Smith } 2503a40ed3dSBarry Smith PetscFunctionReturn(0); 2519b94acceSBarry Smith } 2529b94acceSBarry Smith 253a847f771SSatish Balay 2545615d1e5SSatish Balay #undef __FUNC__ 255d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext" 2569b94acceSBarry Smith /*@ 2579b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 2589b94acceSBarry Smith the nonlinear solvers. 2599b94acceSBarry Smith 2609b94acceSBarry Smith Input Parameters: 2619b94acceSBarry Smith . snes - the SNES context 2629b94acceSBarry Smith . usrP - optional user context 2639b94acceSBarry Smith 2649b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 2659b94acceSBarry Smith 2669b94acceSBarry Smith .seealso: SNESGetApplicationContext() 2679b94acceSBarry Smith @*/ 2689b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 2699b94acceSBarry Smith { 2703a40ed3dSBarry Smith PetscFunctionBegin; 27177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2729b94acceSBarry Smith snes->user = usrP; 2733a40ed3dSBarry Smith PetscFunctionReturn(0); 2749b94acceSBarry Smith } 27574679c65SBarry Smith 2765615d1e5SSatish Balay #undef __FUNC__ 277d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext" 2789b94acceSBarry Smith /*@C 2799b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 2809b94acceSBarry Smith nonlinear solvers. 2819b94acceSBarry Smith 2829b94acceSBarry Smith Input Parameter: 2839b94acceSBarry Smith . snes - SNES context 2849b94acceSBarry Smith 2859b94acceSBarry Smith Output Parameter: 2869b94acceSBarry Smith . usrP - user context 2879b94acceSBarry Smith 2889b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 2899b94acceSBarry Smith 2909b94acceSBarry Smith .seealso: SNESSetApplicationContext() 2919b94acceSBarry Smith @*/ 2929b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 2939b94acceSBarry Smith { 2943a40ed3dSBarry Smith PetscFunctionBegin; 29577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2969b94acceSBarry Smith *usrP = snes->user; 2973a40ed3dSBarry Smith PetscFunctionReturn(0); 2989b94acceSBarry Smith } 29974679c65SBarry Smith 3005615d1e5SSatish Balay #undef __FUNC__ 301d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber" 3029b94acceSBarry Smith /*@ 3039b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3049b94acceSBarry Smith nonlinear solver. 3059b94acceSBarry Smith 3069b94acceSBarry Smith Input Parameter: 3079b94acceSBarry Smith . snes - SNES context 3089b94acceSBarry Smith 3099b94acceSBarry Smith Output Parameter: 3109b94acceSBarry Smith . iter - iteration number 3119b94acceSBarry Smith 3129b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3139b94acceSBarry Smith @*/ 3149b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3159b94acceSBarry Smith { 3163a40ed3dSBarry Smith PetscFunctionBegin; 31777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 31874679c65SBarry Smith PetscValidIntPointer(iter); 3199b94acceSBarry Smith *iter = snes->iter; 3203a40ed3dSBarry Smith PetscFunctionReturn(0); 3219b94acceSBarry Smith } 32274679c65SBarry Smith 3235615d1e5SSatish Balay #undef __FUNC__ 3245615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 3259b94acceSBarry Smith /*@ 3269b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3279b94acceSBarry Smith with SNESSSetFunction(). 3289b94acceSBarry Smith 3299b94acceSBarry Smith Input Parameter: 3309b94acceSBarry Smith . snes - SNES context 3319b94acceSBarry Smith 3329b94acceSBarry Smith Output Parameter: 3339b94acceSBarry Smith . fnorm - 2-norm of function 3349b94acceSBarry Smith 3359b94acceSBarry Smith Note: 3369b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 337a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 338a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 3399b94acceSBarry Smith 3409b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 341a86d99e1SLois Curfman McInnes 342a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 3439b94acceSBarry Smith @*/ 3449b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 3459b94acceSBarry Smith { 3463a40ed3dSBarry Smith PetscFunctionBegin; 34777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 34874679c65SBarry Smith PetscValidScalarPointer(fnorm); 34974679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 350d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 35174679c65SBarry Smith } 3529b94acceSBarry Smith *fnorm = snes->norm; 3533a40ed3dSBarry Smith PetscFunctionReturn(0); 3549b94acceSBarry Smith } 35574679c65SBarry Smith 3565615d1e5SSatish Balay #undef __FUNC__ 3575615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 3589b94acceSBarry Smith /*@ 3599b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 3609b94acceSBarry Smith with SNESSSetGradient(). 3619b94acceSBarry Smith 3629b94acceSBarry Smith Input Parameter: 3639b94acceSBarry Smith . snes - SNES context 3649b94acceSBarry Smith 3659b94acceSBarry Smith Output Parameter: 3669b94acceSBarry Smith . fnorm - 2-norm of gradient 3679b94acceSBarry Smith 3689b94acceSBarry Smith Note: 3699b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 370a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 371a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 3729b94acceSBarry Smith 3739b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 374a86d99e1SLois Curfman McInnes 375a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 3769b94acceSBarry Smith @*/ 3779b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 3789b94acceSBarry Smith { 3793a40ed3dSBarry Smith PetscFunctionBegin; 38077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 38174679c65SBarry Smith PetscValidScalarPointer(gnorm); 38274679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 383d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 38474679c65SBarry Smith } 3859b94acceSBarry Smith *gnorm = snes->norm; 3863a40ed3dSBarry Smith PetscFunctionReturn(0); 3879b94acceSBarry Smith } 38874679c65SBarry Smith 3895615d1e5SSatish Balay #undef __FUNC__ 390d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 3919b94acceSBarry Smith /*@ 3929b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 3939b94acceSBarry Smith attempted by the nonlinear solver. 3949b94acceSBarry Smith 3959b94acceSBarry Smith Input Parameter: 3969b94acceSBarry Smith . snes - SNES context 3979b94acceSBarry Smith 3989b94acceSBarry Smith Output Parameter: 3999b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4009b94acceSBarry Smith 401c96a6f78SLois Curfman McInnes Notes: 402c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 403c96a6f78SLois Curfman McInnes 4049b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4059b94acceSBarry Smith @*/ 4069b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4079b94acceSBarry Smith { 4083a40ed3dSBarry Smith PetscFunctionBegin; 40977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 41074679c65SBarry Smith PetscValidIntPointer(nfails); 4119b94acceSBarry Smith *nfails = snes->nfailures; 4123a40ed3dSBarry Smith PetscFunctionReturn(0); 4139b94acceSBarry Smith } 414a847f771SSatish Balay 4155615d1e5SSatish Balay #undef __FUNC__ 416d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" 417c96a6f78SLois Curfman McInnes /*@ 418c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 419c96a6f78SLois Curfman McInnes used by the nonlinear solver. 420c96a6f78SLois Curfman McInnes 421c96a6f78SLois Curfman McInnes Input Parameter: 422c96a6f78SLois Curfman McInnes . snes - SNES context 423c96a6f78SLois Curfman McInnes 424c96a6f78SLois Curfman McInnes Output Parameter: 425c96a6f78SLois Curfman McInnes . lits - number of linear iterations 426c96a6f78SLois Curfman McInnes 427c96a6f78SLois Curfman McInnes Notes: 428c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 429c96a6f78SLois Curfman McInnes 430c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 431c96a6f78SLois Curfman McInnes @*/ 43286bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 433c96a6f78SLois Curfman McInnes { 4343a40ed3dSBarry Smith PetscFunctionBegin; 435c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 436c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 437c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 4383a40ed3dSBarry Smith PetscFunctionReturn(0); 439c96a6f78SLois Curfman McInnes } 440c96a6f78SLois Curfman McInnes 441c96a6f78SLois Curfman McInnes #undef __FUNC__ 442d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES" 4439b94acceSBarry Smith /*@C 4449b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 4459b94acceSBarry Smith 4469b94acceSBarry Smith Input Parameter: 4479b94acceSBarry Smith . snes - the SNES context 4489b94acceSBarry Smith 4499b94acceSBarry Smith Output Parameter: 4509b94acceSBarry Smith . sles - the SLES context 4519b94acceSBarry Smith 4529b94acceSBarry Smith Notes: 4539b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 4549b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 4559b94acceSBarry Smith KSP and PC contexts as well. 4569b94acceSBarry Smith 4579b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 4589b94acceSBarry Smith 4599b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 4609b94acceSBarry Smith @*/ 4619b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 4629b94acceSBarry Smith { 4633a40ed3dSBarry Smith PetscFunctionBegin; 46477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 4659b94acceSBarry Smith *sles = snes->sles; 4663a40ed3dSBarry Smith PetscFunctionReturn(0); 4679b94acceSBarry Smith } 4689b94acceSBarry Smith /* -----------------------------------------------------------*/ 4695615d1e5SSatish Balay #undef __FUNC__ 4705615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 4719b94acceSBarry Smith /*@C 4729b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 4739b94acceSBarry Smith 4749b94acceSBarry Smith Input Parameter: 4759b94acceSBarry Smith . comm - MPI communicator 4769b94acceSBarry Smith . type - type of method, one of 4779b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS 4789b94acceSBarry Smith $ (for systems of nonlinear equations) 4799b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION 4809b94acceSBarry Smith $ (for unconstrained minimization) 4819b94acceSBarry Smith 4829b94acceSBarry Smith Output Parameter: 4839b94acceSBarry Smith . outsnes - the new SNES context 4849b94acceSBarry Smith 485c1f60f51SBarry Smith Options Database Key: 48619bd6747SLois Curfman McInnes $ -snes_mf - use default matrix-free Jacobian-vector products, 48719bd6747SLois Curfman McInnes $ and no preconditioning matrix 48819bd6747SLois Curfman McInnes $ -snes_mf_operator - use default matrix-free Jacobian-vector 48919bd6747SLois Curfman McInnes $ products, and a user-provided preconditioning matrix 49019bd6747SLois Curfman McInnes $ as set by SNESSetJacobian() 49119bd6747SLois Curfman McInnes $ -snes_fd - use (slow!) finite differences to compute Jacobian 492c1f60f51SBarry Smith 4939b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 4949b94acceSBarry Smith 49563a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 4969b94acceSBarry Smith @*/ 4974b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 4989b94acceSBarry Smith { 4999b94acceSBarry Smith int ierr; 5009b94acceSBarry Smith SNES snes; 5019b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 50237fcc0dbSBarry Smith 5033a40ed3dSBarry Smith PetscFunctionBegin; 5049b94acceSBarry Smith *outsnes = 0; 505d64ed03dSBarry Smith if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){ 506d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 507d64ed03dSBarry Smith } 508d4bb536fSBarry Smith PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm,SNESDestroy,SNESView); 5099b94acceSBarry Smith PLogObjectCreate(snes); 5109b94acceSBarry Smith snes->max_its = 50; 5119750a799SBarry Smith snes->max_funcs = 10000; 5129b94acceSBarry Smith snes->norm = 0.0; 5135a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 514b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 515b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5169b94acceSBarry Smith snes->atol = 1.e-10; 517d64ed03dSBarry Smith } else { 518b4874afaSBarry Smith snes->rtol = 1.e-8; 519b4874afaSBarry Smith snes->ttol = 0.0; 520b4874afaSBarry Smith snes->atol = 1.e-50; 521b4874afaSBarry Smith } 5229b94acceSBarry Smith snes->xtol = 1.e-8; 523b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5249b94acceSBarry Smith snes->nfuncs = 0; 5259b94acceSBarry Smith snes->nfailures = 0; 5267a00f4a9SLois Curfman McInnes snes->linear_its = 0; 527639f9d9dSBarry Smith snes->numbermonitors = 0; 5289b94acceSBarry Smith snes->data = 0; 5299b94acceSBarry Smith snes->view = 0; 5309b94acceSBarry Smith snes->computeumfunction = 0; 5319b94acceSBarry Smith snes->umfunP = 0; 5329b94acceSBarry Smith snes->fc = 0; 5339b94acceSBarry Smith snes->deltatol = 1.e-12; 5349b94acceSBarry Smith snes->fmin = -1.e30; 5359b94acceSBarry Smith snes->method_class = type; 5369b94acceSBarry Smith snes->set_method_called = 0; 537a703fe33SLois Curfman McInnes snes->setup_called = 0; 5389b94acceSBarry Smith snes->ksp_ewconv = 0; 5397d1a2b2bSBarry Smith snes->type = -1; 5406f24a144SLois Curfman McInnes snes->xmonitor = 0; 5416f24a144SLois Curfman McInnes snes->vwork = 0; 5426f24a144SLois Curfman McInnes snes->nwork = 0; 543c9005455SLois Curfman McInnes snes->conv_hist_size = 0; 544c9005455SLois Curfman McInnes snes->conv_act_size = 0; 545c9005455SLois Curfman McInnes snes->conv_hist = 0; 5469b94acceSBarry Smith 5479b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 5480452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 549eed86810SBarry Smith PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 5509b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 5519b94acceSBarry Smith kctx->version = 2; 5529b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 5539b94acceSBarry Smith this was too large for some test cases */ 5549b94acceSBarry Smith kctx->rtol_last = 0; 5559b94acceSBarry Smith kctx->rtol_max = .9; 5569b94acceSBarry Smith kctx->gamma = 1.0; 5579b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 5589b94acceSBarry Smith kctx->alpha = kctx->alpha2; 5599b94acceSBarry Smith kctx->threshold = .1; 5609b94acceSBarry Smith kctx->lresid_last = 0; 5619b94acceSBarry Smith kctx->norm_last = 0; 5629b94acceSBarry Smith 5639b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 5649b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 5655334005bSBarry Smith 5669b94acceSBarry Smith *outsnes = snes; 5673a40ed3dSBarry Smith PetscFunctionReturn(0); 5689b94acceSBarry Smith } 5699b94acceSBarry Smith 5709b94acceSBarry Smith /* --------------------------------------------------------------- */ 5715615d1e5SSatish Balay #undef __FUNC__ 572d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction" 5739b94acceSBarry Smith /*@C 5749b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 5759b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 5769b94acceSBarry Smith equations. 5779b94acceSBarry Smith 5789b94acceSBarry Smith Input Parameters: 5799b94acceSBarry Smith . snes - the SNES context 5809b94acceSBarry Smith . func - function evaluation routine 5819b94acceSBarry Smith . r - vector to store function value 5822cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 5832cd2dfdcSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 5849b94acceSBarry Smith 5859b94acceSBarry Smith Calling sequence of func: 586313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Vec f,void *ctx); 5879b94acceSBarry Smith 5889b94acceSBarry Smith . x - input vector 589313e4042SLois Curfman McInnes . f - function vector 5902cd2dfdcSLois Curfman McInnes . ctx - optional user-defined function context 5919b94acceSBarry Smith 5929b94acceSBarry Smith Notes: 5939b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 5949b94acceSBarry Smith $ f'(x) x = -f(x), 5959b94acceSBarry Smith $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 5969b94acceSBarry Smith 5979b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 5989b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 5999b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 6009b94acceSBarry Smith 6019b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6029b94acceSBarry Smith 603a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6049b94acceSBarry Smith @*/ 6055334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 6069b94acceSBarry Smith { 6073a40ed3dSBarry Smith PetscFunctionBegin; 60877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 609*a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 610*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 611*a8c6a408SBarry Smith } 6129b94acceSBarry Smith snes->computefunction = func; 6139b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6149b94acceSBarry Smith snes->funP = ctx; 6153a40ed3dSBarry Smith PetscFunctionReturn(0); 6169b94acceSBarry Smith } 6179b94acceSBarry Smith 6185615d1e5SSatish Balay #undef __FUNC__ 6195615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 6209b94acceSBarry Smith /*@ 6219b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6229b94acceSBarry Smith SNESSetFunction(). 6239b94acceSBarry Smith 6249b94acceSBarry Smith Input Parameters: 6259b94acceSBarry Smith . snes - the SNES context 6269b94acceSBarry Smith . x - input vector 6279b94acceSBarry Smith 6289b94acceSBarry Smith Output Parameter: 6293638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6309b94acceSBarry Smith 6311bffabb2SLois Curfman McInnes Notes: 6329b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6339b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6349b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6359b94acceSBarry Smith 6369b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6379b94acceSBarry Smith 638a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6399b94acceSBarry Smith @*/ 6409b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 6419b94acceSBarry Smith { 6429b94acceSBarry Smith int ierr; 6439b94acceSBarry Smith 6443a40ed3dSBarry Smith PetscFunctionBegin; 645d4bb536fSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 646*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 647d4bb536fSBarry Smith } 6489b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 649d64ed03dSBarry Smith PetscStackPush("SNES user function"); 6509b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 651d64ed03dSBarry Smith PetscStackPop; 652ae3c334cSLois Curfman McInnes snes->nfuncs++; 6539b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 6543a40ed3dSBarry Smith PetscFunctionReturn(0); 6559b94acceSBarry Smith } 6569b94acceSBarry Smith 6575615d1e5SSatish Balay #undef __FUNC__ 658d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" 6599b94acceSBarry Smith /*@C 6609b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 6619b94acceSBarry Smith unconstrained minimization. 6629b94acceSBarry Smith 6639b94acceSBarry Smith Input Parameters: 6649b94acceSBarry Smith . snes - the SNES context 6659b94acceSBarry Smith . func - function evaluation routine 666044dda88SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 667044dda88SLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6689b94acceSBarry Smith 6699b94acceSBarry Smith Calling sequence of func: 6709b94acceSBarry Smith . func (SNES snes,Vec x,double *f,void *ctx); 6719b94acceSBarry Smith 6729b94acceSBarry Smith . x - input vector 6739b94acceSBarry Smith . f - function 674044dda88SLois Curfman McInnes . ctx - [optional] user-defined function context 6759b94acceSBarry Smith 6769b94acceSBarry Smith Notes: 6779b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 6789b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 6799b94acceSBarry Smith SNESSetFunction(). 6809b94acceSBarry Smith 6819b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 6829b94acceSBarry Smith 683a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 684a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 6859b94acceSBarry Smith @*/ 6869b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 6879b94acceSBarry Smith void *ctx) 6889b94acceSBarry Smith { 6893a40ed3dSBarry Smith PetscFunctionBegin; 69077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 691*a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 692*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 693*a8c6a408SBarry Smith } 6949b94acceSBarry Smith snes->computeumfunction = func; 6959b94acceSBarry Smith snes->umfunP = ctx; 6963a40ed3dSBarry Smith PetscFunctionReturn(0); 6979b94acceSBarry Smith } 6989b94acceSBarry Smith 6995615d1e5SSatish Balay #undef __FUNC__ 7005615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 7019b94acceSBarry Smith /*@ 7029b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 7039b94acceSBarry Smith set with SNESSetMinimizationFunction(). 7049b94acceSBarry Smith 7059b94acceSBarry Smith Input Parameters: 7069b94acceSBarry Smith . snes - the SNES context 7079b94acceSBarry Smith . x - input vector 7089b94acceSBarry Smith 7099b94acceSBarry Smith Output Parameter: 7109b94acceSBarry Smith . y - function value 7119b94acceSBarry Smith 7129b94acceSBarry Smith Notes: 7139b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 7149b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7159b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 716a86d99e1SLois Curfman McInnes 717a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 718a86d99e1SLois Curfman McInnes 719a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 720a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 7219b94acceSBarry Smith @*/ 7229b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 7239b94acceSBarry Smith { 7249b94acceSBarry Smith int ierr; 7253a40ed3dSBarry Smith 7263a40ed3dSBarry Smith PetscFunctionBegin; 727*a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 728*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 729*a8c6a408SBarry Smith } 7309b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 731d64ed03dSBarry Smith PetscStackPush("SNES user minimzation function"); 7329b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 733d64ed03dSBarry Smith PetscStackPop; 734ae3c334cSLois Curfman McInnes snes->nfuncs++; 7359b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7363a40ed3dSBarry Smith PetscFunctionReturn(0); 7379b94acceSBarry Smith } 7389b94acceSBarry Smith 7395615d1e5SSatish Balay #undef __FUNC__ 740d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient" 7419b94acceSBarry Smith /*@C 7429b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 7439b94acceSBarry Smith vector for use by the SNES routines. 7449b94acceSBarry Smith 7459b94acceSBarry Smith Input Parameters: 7469b94acceSBarry Smith . snes - the SNES context 7479b94acceSBarry Smith . func - function evaluation routine 748044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 749044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 7509b94acceSBarry Smith . r - vector to store gradient value 7519b94acceSBarry Smith 7529b94acceSBarry Smith Calling sequence of func: 7539b94acceSBarry Smith . func (SNES, Vec x, Vec g, void *ctx); 7549b94acceSBarry Smith 7559b94acceSBarry Smith . x - input vector 7569b94acceSBarry Smith . g - gradient vector 757044dda88SLois Curfman McInnes . ctx - optional user-defined gradient context 7589b94acceSBarry Smith 7599b94acceSBarry Smith Notes: 7609b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7619b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7629b94acceSBarry Smith SNESSetFunction(). 7639b94acceSBarry Smith 7649b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 7659b94acceSBarry Smith 766a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 767a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 7689b94acceSBarry Smith @*/ 76974679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 7709b94acceSBarry Smith { 7713a40ed3dSBarry Smith PetscFunctionBegin; 77277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 773*a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 774*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 775*a8c6a408SBarry Smith } 7769b94acceSBarry Smith snes->computefunction = func; 7779b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7789b94acceSBarry Smith snes->funP = ctx; 7793a40ed3dSBarry Smith PetscFunctionReturn(0); 7809b94acceSBarry Smith } 7819b94acceSBarry Smith 7825615d1e5SSatish Balay #undef __FUNC__ 7835615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 7849b94acceSBarry Smith /*@ 785a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 786a86d99e1SLois Curfman McInnes SNESSetGradient(). 7879b94acceSBarry Smith 7889b94acceSBarry Smith Input Parameters: 7899b94acceSBarry Smith . snes - the SNES context 7909b94acceSBarry Smith . x - input vector 7919b94acceSBarry Smith 7929b94acceSBarry Smith Output Parameter: 7939b94acceSBarry Smith . y - gradient vector 7949b94acceSBarry Smith 7959b94acceSBarry Smith Notes: 7969b94acceSBarry Smith SNESComputeGradient() is valid only for 7979b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7989b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 799a86d99e1SLois Curfman McInnes 800a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 801a86d99e1SLois Curfman McInnes 802a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 803a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 8049b94acceSBarry Smith @*/ 8059b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 8069b94acceSBarry Smith { 8079b94acceSBarry Smith int ierr; 8083a40ed3dSBarry Smith 8093a40ed3dSBarry Smith PetscFunctionBegin; 8103a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 811*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8123a40ed3dSBarry Smith } 8139b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 814d64ed03dSBarry Smith PetscStackPush("SNES user gradient function"); 8159b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 816d64ed03dSBarry Smith PetscStackPop; 8179b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 8183a40ed3dSBarry Smith PetscFunctionReturn(0); 8199b94acceSBarry Smith } 8209b94acceSBarry Smith 8215615d1e5SSatish Balay #undef __FUNC__ 8225615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 82362fef451SLois Curfman McInnes /*@ 82462fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 82562fef451SLois Curfman McInnes set with SNESSetJacobian(). 82662fef451SLois Curfman McInnes 82762fef451SLois Curfman McInnes Input Parameters: 82862fef451SLois Curfman McInnes . snes - the SNES context 82962fef451SLois Curfman McInnes . x - input vector 83062fef451SLois Curfman McInnes 83162fef451SLois Curfman McInnes Output Parameters: 83262fef451SLois Curfman McInnes . A - Jacobian matrix 83362fef451SLois Curfman McInnes . B - optional preconditioning matrix 83462fef451SLois Curfman McInnes . flag - flag indicating matrix structure 83562fef451SLois Curfman McInnes 83662fef451SLois Curfman McInnes Notes: 83762fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 83862fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 83962fef451SLois Curfman McInnes 840dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 841dc5a77f8SLois Curfman McInnes flag parameter. 84262fef451SLois Curfman McInnes 84362fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 84462fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 845005c665bSBarry Smith methods is SNESComputeHessian(). 84662fef451SLois Curfman McInnes 84762fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 84862fef451SLois Curfman McInnes 84962fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 85062fef451SLois Curfman McInnes @*/ 8519b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8529b94acceSBarry Smith { 8539b94acceSBarry Smith int ierr; 8543a40ed3dSBarry Smith 8553a40ed3dSBarry Smith PetscFunctionBegin; 8563a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 857*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 8583a40ed3dSBarry Smith } 8593a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 8609b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 861c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 862d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 8639b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 864d64ed03dSBarry Smith PetscStackPop; 8659b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 8666d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 86777c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 86877c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 8693a40ed3dSBarry Smith PetscFunctionReturn(0); 8709b94acceSBarry Smith } 8719b94acceSBarry Smith 8725615d1e5SSatish Balay #undef __FUNC__ 8735615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 87462fef451SLois Curfman McInnes /*@ 87562fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 87662fef451SLois Curfman McInnes set with SNESSetHessian(). 87762fef451SLois Curfman McInnes 87862fef451SLois Curfman McInnes Input Parameters: 87962fef451SLois Curfman McInnes . snes - the SNES context 88062fef451SLois Curfman McInnes . x - input vector 88162fef451SLois Curfman McInnes 88262fef451SLois Curfman McInnes Output Parameters: 88362fef451SLois Curfman McInnes . A - Hessian matrix 88462fef451SLois Curfman McInnes . B - optional preconditioning matrix 88562fef451SLois Curfman McInnes . flag - flag indicating matrix structure 88662fef451SLois Curfman McInnes 88762fef451SLois Curfman McInnes Notes: 88862fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 88962fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 89062fef451SLois Curfman McInnes 891dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 892dc5a77f8SLois Curfman McInnes flag parameter. 89362fef451SLois Curfman McInnes 89462fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 89562fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 89662fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 89762fef451SLois Curfman McInnes 89862fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 89962fef451SLois Curfman McInnes 900a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 901a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 90262fef451SLois Curfman McInnes @*/ 90362fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 9049b94acceSBarry Smith { 9059b94acceSBarry Smith int ierr; 9063a40ed3dSBarry Smith 9073a40ed3dSBarry Smith PetscFunctionBegin; 9083a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 909*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 9103a40ed3dSBarry Smith } 9113a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 91262fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 9130f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 914d64ed03dSBarry Smith PetscStackPush("SNES user Hessian function"); 91562fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 916d64ed03dSBarry Smith PetscStackPop; 91762fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 9180f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 91977c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 92077c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9213a40ed3dSBarry Smith PetscFunctionReturn(0); 9229b94acceSBarry Smith } 9239b94acceSBarry Smith 9245615d1e5SSatish Balay #undef __FUNC__ 925d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian" 9269b94acceSBarry Smith /*@C 9279b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 928044dda88SLois Curfman McInnes location to store the matrix. 9299b94acceSBarry Smith 9309b94acceSBarry Smith Input Parameters: 9319b94acceSBarry Smith . snes - the SNES context 9329b94acceSBarry Smith . A - Jacobian matrix 9339b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 9349b94acceSBarry Smith . func - Jacobian evaluation routine 9352cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 9362cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 9379b94acceSBarry Smith 9389b94acceSBarry Smith Calling sequence of func: 939313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9409b94acceSBarry Smith 9419b94acceSBarry Smith . x - input vector 9429b94acceSBarry Smith . A - Jacobian matrix 9439b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 944ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 945ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 9462cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Jacobian context 9479b94acceSBarry Smith 9489b94acceSBarry Smith Notes: 949dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 9502cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 951ac21db08SLois Curfman McInnes 952ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9539b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9549b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9559b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9569b94acceSBarry Smith throughout the global iterations. 9579b94acceSBarry Smith 9589b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9599b94acceSBarry Smith 960ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 9619b94acceSBarry Smith @*/ 9629b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 9639b94acceSBarry Smith MatStructure*,void*),void *ctx) 9649b94acceSBarry Smith { 9653a40ed3dSBarry Smith PetscFunctionBegin; 96677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 967*a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 968*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 969*a8c6a408SBarry Smith } 9709b94acceSBarry Smith snes->computejacobian = func; 9719b94acceSBarry Smith snes->jacP = ctx; 9729b94acceSBarry Smith snes->jacobian = A; 9739b94acceSBarry Smith snes->jacobian_pre = B; 9743a40ed3dSBarry Smith PetscFunctionReturn(0); 9759b94acceSBarry Smith } 97662fef451SLois Curfman McInnes 9775615d1e5SSatish Balay #undef __FUNC__ 978d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian" 979b4fd4287SBarry Smith /*@ 980b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 981b4fd4287SBarry Smith provided context for evaluating the Jacobian. 982b4fd4287SBarry Smith 983b4fd4287SBarry Smith Input Parameter: 984b4fd4287SBarry Smith . snes - the nonlinear solver context 985b4fd4287SBarry Smith 986b4fd4287SBarry Smith Output Parameters: 987b4fd4287SBarry Smith . A - location to stash Jacobian matrix (or PETSC_NULL) 988b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 989b4fd4287SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 990b4fd4287SBarry Smith 991b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 992b4fd4287SBarry Smith @*/ 993b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 994b4fd4287SBarry Smith { 9953a40ed3dSBarry Smith PetscFunctionBegin; 9963a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 997*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 9983a40ed3dSBarry Smith } 999b4fd4287SBarry Smith if (A) *A = snes->jacobian; 1000b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 1001b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 10023a40ed3dSBarry Smith PetscFunctionReturn(0); 1003b4fd4287SBarry Smith } 1004b4fd4287SBarry Smith 10055615d1e5SSatish Balay #undef __FUNC__ 1006d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian" 10079b94acceSBarry Smith /*@C 10089b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 1009044dda88SLois Curfman McInnes location to store the matrix. 10109b94acceSBarry Smith 10119b94acceSBarry Smith Input Parameters: 10129b94acceSBarry Smith . snes - the SNES context 10139b94acceSBarry Smith . A - Hessian matrix 10149b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 10159b94acceSBarry Smith . func - Jacobian evaluation routine 1016313e4042SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 1017313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 10189b94acceSBarry Smith 10199b94acceSBarry Smith Calling sequence of func: 1020313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 10219b94acceSBarry Smith 10229b94acceSBarry Smith . x - input vector 10239b94acceSBarry Smith . A - Hessian matrix 10249b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1025ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1026ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 10272cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Hessian context 10289b94acceSBarry Smith 10299b94acceSBarry Smith Notes: 1030dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10312cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1032ac21db08SLois Curfman McInnes 10339b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 10349b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 10359b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10369b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10379b94acceSBarry Smith throughout the global iterations. 10389b94acceSBarry Smith 10399b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 10409b94acceSBarry Smith 1041ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 10429b94acceSBarry Smith @*/ 10439b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 10449b94acceSBarry Smith MatStructure*,void*),void *ctx) 10459b94acceSBarry Smith { 10463a40ed3dSBarry Smith PetscFunctionBegin; 104777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1048d4bb536fSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1049*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1050d4bb536fSBarry Smith } 10519b94acceSBarry Smith snes->computejacobian = func; 10529b94acceSBarry Smith snes->jacP = ctx; 10539b94acceSBarry Smith snes->jacobian = A; 10549b94acceSBarry Smith snes->jacobian_pre = B; 10553a40ed3dSBarry Smith PetscFunctionReturn(0); 10569b94acceSBarry Smith } 10579b94acceSBarry Smith 10585615d1e5SSatish Balay #undef __FUNC__ 1059d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian" 106062fef451SLois Curfman McInnes /*@ 106162fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 106262fef451SLois Curfman McInnes provided context for evaluating the Hessian. 106362fef451SLois Curfman McInnes 106462fef451SLois Curfman McInnes Input Parameter: 106562fef451SLois Curfman McInnes . snes - the nonlinear solver context 106662fef451SLois Curfman McInnes 106762fef451SLois Curfman McInnes Output Parameters: 106862fef451SLois Curfman McInnes . A - location to stash Hessian matrix (or PETSC_NULL) 106962fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 107062fef451SLois Curfman McInnes . ctx - location to stash Hessian ctx (or PETSC_NULL) 107162fef451SLois Curfman McInnes 107262fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 107362fef451SLois Curfman McInnes @*/ 107462fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 107562fef451SLois Curfman McInnes { 10763a40ed3dSBarry Smith PetscFunctionBegin; 10773a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){ 1078*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 10793a40ed3dSBarry Smith } 108062fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 108162fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 108262fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 10833a40ed3dSBarry Smith PetscFunctionReturn(0); 108462fef451SLois Curfman McInnes } 108562fef451SLois Curfman McInnes 10869b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 10879b94acceSBarry Smith 10885615d1e5SSatish Balay #undef __FUNC__ 10895615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 10909b94acceSBarry Smith /*@ 10919b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1092272ac6f2SLois Curfman McInnes of a nonlinear solver. 10939b94acceSBarry Smith 10949b94acceSBarry Smith Input Parameter: 10959b94acceSBarry Smith . snes - the SNES context 10968ddd3da0SLois Curfman McInnes . x - the solution vector 10979b94acceSBarry Smith 1098272ac6f2SLois Curfman McInnes Notes: 1099272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1100272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1101272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1102272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1103272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1104272ac6f2SLois Curfman McInnes 11059b94acceSBarry Smith .keywords: SNES, nonlinear, setup 11069b94acceSBarry Smith 11079b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 11089b94acceSBarry Smith @*/ 11098ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 11109b94acceSBarry Smith { 1111272ac6f2SLois Curfman McInnes int ierr, flg; 11123a40ed3dSBarry Smith 11133a40ed3dSBarry Smith PetscFunctionBegin; 111477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 111577c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 11168ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 11179b94acceSBarry Smith 1118c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1119c1f60f51SBarry Smith /* 1120c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1121dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1122c1f60f51SBarry Smith */ 1123c1f60f51SBarry Smith if (flg) { 1124c1f60f51SBarry Smith Mat J; 1125c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1126c1f60f51SBarry Smith PLogObjectParent(snes,J); 1127c1f60f51SBarry Smith snes->mfshell = J; 1128c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1129c1f60f51SBarry Smith snes->jacobian = J; 1130c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1131d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1132c1f60f51SBarry Smith snes->jacobian = J; 1133c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1134d4bb536fSBarry Smith } else { 1135*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option"); 1136d4bb536fSBarry Smith } 1137c1f60f51SBarry Smith } 1138272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1139c1f60f51SBarry Smith /* 1140dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1141c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1142c1f60f51SBarry Smith */ 114331615d04SLois Curfman McInnes if (flg) { 1144272ac6f2SLois Curfman McInnes Mat J; 1145272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1146272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1147272ac6f2SLois Curfman McInnes snes->mfshell = J; 114883e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 114983e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 115094a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1151d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 115283e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 115394a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1154d4bb536fSBarry Smith } else { 1155*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option"); 1156d4bb536fSBarry Smith } 1157272ac6f2SLois Curfman McInnes } 11589b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1159*a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1160*a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1161*a8c6a408SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first"); 1162*a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 1163*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector"); 1164*a8c6a408SBarry Smith } 1165a703fe33SLois Curfman McInnes 1166a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 116740191667SLois Curfman McInnes if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1168a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1169a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1170a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1171a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1172a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1173a703fe33SLois Curfman McInnes } 11749b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1175*a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1176*a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1177*a8c6a408SBarry Smith if (!snes->computeumfunction) { 1178*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first"); 1179*a8c6a408SBarry Smith } 1180*a8c6a408SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first"); 1181d4bb536fSBarry Smith } else { 1182*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class"); 1183d4bb536fSBarry Smith } 1184a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1185a703fe33SLois Curfman McInnes snes->setup_called = 1; 11863a40ed3dSBarry Smith PetscFunctionReturn(0); 11879b94acceSBarry Smith } 11889b94acceSBarry Smith 11895615d1e5SSatish Balay #undef __FUNC__ 1190d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy" 11919b94acceSBarry Smith /*@C 11929b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 11939b94acceSBarry Smith with SNESCreate(). 11949b94acceSBarry Smith 11959b94acceSBarry Smith Input Parameter: 11969b94acceSBarry Smith . snes - the SNES context 11979b94acceSBarry Smith 11989b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 11999b94acceSBarry Smith 120063a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 12019b94acceSBarry Smith @*/ 12029b94acceSBarry Smith int SNESDestroy(SNES snes) 12039b94acceSBarry Smith { 12049b94acceSBarry Smith int ierr; 12053a40ed3dSBarry Smith 12063a40ed3dSBarry Smith PetscFunctionBegin; 120777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12083a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1209d4bb536fSBarry Smith 12109750a799SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);} 12110452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 1212522c5e43SBarry Smith if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);} 12139b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 1214522c5e43SBarry Smith if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);} 1215522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 12169b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 12170452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 12183a40ed3dSBarry Smith PetscFunctionReturn(0); 12199b94acceSBarry Smith } 12209b94acceSBarry Smith 12219b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 12229b94acceSBarry Smith 12235615d1e5SSatish Balay #undef __FUNC__ 12245615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 12259b94acceSBarry Smith /*@ 1226d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 12279b94acceSBarry Smith 12289b94acceSBarry Smith Input Parameters: 12299b94acceSBarry Smith . snes - the SNES context 123033174efeSLois Curfman McInnes . atol - absolute convergence tolerance 123133174efeSLois Curfman McInnes . rtol - relative convergence tolerance 123233174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 123333174efeSLois Curfman McInnes of the change in the solution between steps 123433174efeSLois Curfman McInnes . maxit - maximum number of iterations 123533174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 12369b94acceSBarry Smith 123733174efeSLois Curfman McInnes Options Database Keys: 123833174efeSLois Curfman McInnes $ -snes_atol <atol> 123933174efeSLois Curfman McInnes $ -snes_rtol <rtol> 124033174efeSLois Curfman McInnes $ -snes_stol <stol> 124133174efeSLois Curfman McInnes $ -snes_max_it <maxit> 124233174efeSLois Curfman McInnes $ -snes_max_funcs <maxf> 12439b94acceSBarry Smith 1244d7a720efSLois Curfman McInnes Notes: 12459b94acceSBarry Smith The default maximum number of iterations is 50. 12469b94acceSBarry Smith The default maximum number of function evaluations is 1000. 12479b94acceSBarry Smith 124833174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 12499b94acceSBarry Smith 1250d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 12519b94acceSBarry Smith @*/ 1252d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 12539b94acceSBarry Smith { 12543a40ed3dSBarry Smith PetscFunctionBegin; 125577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1256d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1257d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1258d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1259d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1260d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 12613a40ed3dSBarry Smith PetscFunctionReturn(0); 12629b94acceSBarry Smith } 12639b94acceSBarry Smith 12645615d1e5SSatish Balay #undef __FUNC__ 12655615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 12669b94acceSBarry Smith /*@ 126733174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 126833174efeSLois Curfman McInnes 126933174efeSLois Curfman McInnes Input Parameters: 127033174efeSLois Curfman McInnes . snes - the SNES context 127133174efeSLois Curfman McInnes . atol - absolute convergence tolerance 127233174efeSLois Curfman McInnes . rtol - relative convergence tolerance 127333174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 127433174efeSLois Curfman McInnes of the change in the solution between steps 127533174efeSLois Curfman McInnes . maxit - maximum number of iterations 127633174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 127733174efeSLois Curfman McInnes 127833174efeSLois Curfman McInnes Notes: 127933174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 128033174efeSLois Curfman McInnes 128133174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 128233174efeSLois Curfman McInnes 128333174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 128433174efeSLois Curfman McInnes @*/ 128533174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 128633174efeSLois Curfman McInnes { 12873a40ed3dSBarry Smith PetscFunctionBegin; 128833174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 128933174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 129033174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 129133174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 129233174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 129333174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 12943a40ed3dSBarry Smith PetscFunctionReturn(0); 129533174efeSLois Curfman McInnes } 129633174efeSLois Curfman McInnes 12975615d1e5SSatish Balay #undef __FUNC__ 12985615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 129933174efeSLois Curfman McInnes /*@ 13009b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 13019b94acceSBarry Smith 13029b94acceSBarry Smith Input Parameters: 13039b94acceSBarry Smith . snes - the SNES context 13049b94acceSBarry Smith . tol - tolerance 13059b94acceSBarry Smith 13069b94acceSBarry Smith Options Database Key: 1307d7a720efSLois Curfman McInnes $ -snes_trtol <tol> 13089b94acceSBarry Smith 13099b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 13109b94acceSBarry Smith 1311d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 13129b94acceSBarry Smith @*/ 13139b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 13149b94acceSBarry Smith { 13153a40ed3dSBarry Smith PetscFunctionBegin; 131677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13179b94acceSBarry Smith snes->deltatol = tol; 13183a40ed3dSBarry Smith PetscFunctionReturn(0); 13199b94acceSBarry Smith } 13209b94acceSBarry Smith 13215615d1e5SSatish Balay #undef __FUNC__ 13225615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 13239b94acceSBarry Smith /*@ 132477c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 13259b94acceSBarry Smith for unconstrained minimization solvers. 13269b94acceSBarry Smith 13279b94acceSBarry Smith Input Parameters: 13289b94acceSBarry Smith . snes - the SNES context 13299b94acceSBarry Smith . ftol - minimum function tolerance 13309b94acceSBarry Smith 13319b94acceSBarry Smith Options Database Key: 1332d7a720efSLois Curfman McInnes $ -snes_fmin <ftol> 13339b94acceSBarry Smith 13349b94acceSBarry Smith Note: 133577c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 13369b94acceSBarry Smith methods only. 13379b94acceSBarry Smith 13389b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 13399b94acceSBarry Smith 1340d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 13419b94acceSBarry Smith @*/ 134277c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 13439b94acceSBarry Smith { 13443a40ed3dSBarry Smith PetscFunctionBegin; 134577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13469b94acceSBarry Smith snes->fmin = ftol; 13473a40ed3dSBarry Smith PetscFunctionReturn(0); 13489b94acceSBarry Smith } 13499b94acceSBarry Smith 13509b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 13519b94acceSBarry Smith 13525615d1e5SSatish Balay #undef __FUNC__ 1353d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor" 13549b94acceSBarry Smith /*@C 13559b94acceSBarry Smith SNESSetMonitor - Sets the function that is to be used at every 13569b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 13579b94acceSBarry Smith progress. 13589b94acceSBarry Smith 13599b94acceSBarry Smith Input Parameters: 13609b94acceSBarry Smith . snes - the SNES context 13619b94acceSBarry Smith . func - monitoring routine 1362044dda88SLois Curfman McInnes . mctx - [optional] user-defined context for private data for the 1363044dda88SLois Curfman McInnes monitor routine (may be PETSC_NULL) 13649b94acceSBarry Smith 13659b94acceSBarry Smith Calling sequence of func: 1366682d7d0cSBarry Smith int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 13679b94acceSBarry Smith 13689b94acceSBarry Smith $ snes - the SNES context 13699b94acceSBarry Smith $ its - iteration number 1370044dda88SLois Curfman McInnes $ mctx - [optional] monitoring context 13719b94acceSBarry Smith $ 13729b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13739b94acceSBarry Smith $ norm - 2-norm function value (may be estimated) 13749b94acceSBarry Smith $ 13759b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13769b94acceSBarry Smith $ norm - 2-norm gradient value (may be estimated) 13779b94acceSBarry Smith 13789665c990SLois Curfman McInnes Options Database Keys: 13799665c990SLois Curfman McInnes $ -snes_monitor : sets SNESDefaultMonitor() 13809665c990SLois Curfman McInnes $ -snes_xmonitor : sets line graph monitor, 13819665c990SLois Curfman McInnes $ uses SNESLGMonitorCreate() 13829665c990SLois Curfman McInnes $ -snes_cancelmonitors : cancels all monitors that have 13839665c990SLois Curfman McInnes $ been hardwired into a code by 13849665c990SLois Curfman McInnes $ calls to SNESSetMonitor(), but 13859665c990SLois Curfman McInnes $ does not cancel those set via 13869665c990SLois Curfman McInnes $ the options database. 13879665c990SLois Curfman McInnes 13889665c990SLois Curfman McInnes 1389639f9d9dSBarry Smith Notes: 13906bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 13916bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 13926bc08f3fSLois Curfman McInnes order in which they were set. 1393639f9d9dSBarry Smith 13949b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 13959b94acceSBarry Smith 13969b94acceSBarry Smith .seealso: SNESDefaultMonitor() 13979b94acceSBarry Smith @*/ 139874679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 13999b94acceSBarry Smith { 14003a40ed3dSBarry Smith PetscFunctionBegin; 1401639f9d9dSBarry Smith if (!func) { 1402639f9d9dSBarry Smith snes->numbermonitors = 0; 14033a40ed3dSBarry Smith PetscFunctionReturn(0); 1404639f9d9dSBarry Smith } 1405639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1406*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set"); 1407639f9d9dSBarry Smith } 1408639f9d9dSBarry Smith 1409639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1410639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 14113a40ed3dSBarry Smith PetscFunctionReturn(0); 14129b94acceSBarry Smith } 14139b94acceSBarry Smith 14145615d1e5SSatish Balay #undef __FUNC__ 1415d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest" 14169b94acceSBarry Smith /*@C 14179b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 14189b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 14199b94acceSBarry Smith 14209b94acceSBarry Smith Input Parameters: 14219b94acceSBarry Smith . snes - the SNES context 14229b94acceSBarry Smith . func - routine to test for convergence 1423044dda88SLois Curfman McInnes . cctx - [optional] context for private data for the convergence routine 1424044dda88SLois Curfman McInnes (may be PETSC_NULL) 14259b94acceSBarry Smith 14269b94acceSBarry Smith Calling sequence of func: 14279b94acceSBarry Smith int func (SNES snes,double xnorm,double gnorm, 14289b94acceSBarry Smith double f,void *cctx) 14299b94acceSBarry Smith 14309b94acceSBarry Smith $ snes - the SNES context 1431044dda88SLois Curfman McInnes $ cctx - [optional] convergence context 14329b94acceSBarry Smith $ xnorm - 2-norm of current iterate 14339b94acceSBarry Smith $ 14349b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 14359b94acceSBarry Smith $ gnorm - 2-norm of current step 14369b94acceSBarry Smith $ f - 2-norm of function 14379b94acceSBarry Smith $ 14389b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 14399b94acceSBarry Smith $ gnorm - 2-norm of current gradient 14409b94acceSBarry Smith $ f - function value 14419b94acceSBarry Smith 14429b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 14439b94acceSBarry Smith 144440191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 144540191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 14469b94acceSBarry Smith @*/ 144774679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 14489b94acceSBarry Smith { 14493a40ed3dSBarry Smith PetscFunctionBegin; 14509b94acceSBarry Smith (snes)->converged = func; 14519b94acceSBarry Smith (snes)->cnvP = cctx; 14523a40ed3dSBarry Smith PetscFunctionReturn(0); 14539b94acceSBarry Smith } 14549b94acceSBarry Smith 14555615d1e5SSatish Balay #undef __FUNC__ 1456d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" 1457c9005455SLois Curfman McInnes /*@ 1458c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1459c9005455SLois Curfman McInnes 1460c9005455SLois Curfman McInnes Input Parameters: 1461c9005455SLois Curfman McInnes . snes - iterative context obtained from SNESCreate() 1462c9005455SLois Curfman McInnes . a - array to hold history 1463c9005455SLois Curfman McInnes . na - size of a 1464c9005455SLois Curfman McInnes 1465c9005455SLois Curfman McInnes Notes: 1466c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1467c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1468c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1469c9005455SLois Curfman McInnes at each step. 1470c9005455SLois Curfman McInnes 1471c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1472c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1473c9005455SLois Curfman McInnes during the section of code that is being timed. 1474c9005455SLois Curfman McInnes 1475c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1476c9005455SLois Curfman McInnes @*/ 1477c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na) 1478c9005455SLois Curfman McInnes { 14793a40ed3dSBarry Smith PetscFunctionBegin; 1480c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1481c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1482c9005455SLois Curfman McInnes snes->conv_hist = a; 1483c9005455SLois Curfman McInnes snes->conv_hist_size = na; 14843a40ed3dSBarry Smith PetscFunctionReturn(0); 1485c9005455SLois Curfman McInnes } 1486c9005455SLois Curfman McInnes 1487c9005455SLois Curfman McInnes #undef __FUNC__ 14885615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 14899b94acceSBarry Smith /* 14909b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 14919b94acceSBarry Smith positive parameter delta. 14929b94acceSBarry Smith 14939b94acceSBarry Smith Input Parameters: 14949b94acceSBarry Smith . snes - the SNES context 14959b94acceSBarry Smith . y - approximate solution of linear system 14969b94acceSBarry Smith . fnorm - 2-norm of current function 14979b94acceSBarry Smith . delta - trust region size 14989b94acceSBarry Smith 14999b94acceSBarry Smith Output Parameters: 15009b94acceSBarry Smith . gpnorm - predicted function norm at the new point, assuming local 15019b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 15029b94acceSBarry Smith region, and exceeds zero otherwise. 15039b94acceSBarry Smith . ynorm - 2-norm of the step 15049b94acceSBarry Smith 15059b94acceSBarry Smith Note: 150640191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 15079b94acceSBarry Smith is set to be the maximum allowable step size. 15089b94acceSBarry Smith 15099b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 15109b94acceSBarry Smith */ 15119b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 15129b94acceSBarry Smith double *gpnorm,double *ynorm) 15139b94acceSBarry Smith { 15149b94acceSBarry Smith double norm; 15159b94acceSBarry Smith Scalar cnorm; 15163a40ed3dSBarry Smith int ierr; 15173a40ed3dSBarry Smith 15183a40ed3dSBarry Smith PetscFunctionBegin; 15193a40ed3dSBarry Smith ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr); 15209b94acceSBarry Smith if (norm > *delta) { 15219b94acceSBarry Smith norm = *delta/norm; 15229b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 15239b94acceSBarry Smith cnorm = norm; 15249b94acceSBarry Smith VecScale( &cnorm, y ); 15259b94acceSBarry Smith *ynorm = *delta; 15269b94acceSBarry Smith } else { 15279b94acceSBarry Smith *gpnorm = 0.0; 15289b94acceSBarry Smith *ynorm = norm; 15299b94acceSBarry Smith } 15303a40ed3dSBarry Smith PetscFunctionReturn(0); 15319b94acceSBarry Smith } 15329b94acceSBarry Smith 15335615d1e5SSatish Balay #undef __FUNC__ 15345615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 15359b94acceSBarry Smith /*@ 15369b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 153763a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 15389b94acceSBarry Smith 15399b94acceSBarry Smith Input Parameter: 15409b94acceSBarry Smith . snes - the SNES context 15418ddd3da0SLois Curfman McInnes . x - the solution vector 15429b94acceSBarry Smith 15439b94acceSBarry Smith Output Parameter: 15449b94acceSBarry Smith its - number of iterations until termination 15459b94acceSBarry Smith 15468ddd3da0SLois Curfman McInnes Note: 15478ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 15488ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 15498ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 15508ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 15518ddd3da0SLois Curfman McInnes 15529b94acceSBarry Smith .keywords: SNES, nonlinear, solve 15539b94acceSBarry Smith 155463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 15559b94acceSBarry Smith @*/ 15568ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 15579b94acceSBarry Smith { 15583c7409f5SSatish Balay int ierr, flg; 1559052efed2SBarry Smith 15603a40ed3dSBarry Smith PetscFunctionBegin; 156177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 156274679c65SBarry Smith PetscValidIntPointer(its); 1563c4fc05e7SBarry Smith if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1564c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 15659b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 1566c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 15679b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 15689b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 15693c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 15706d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 15713a40ed3dSBarry Smith PetscFunctionReturn(0); 15729b94acceSBarry Smith } 15739b94acceSBarry Smith 15749b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 15759b94acceSBarry Smith 15765615d1e5SSatish Balay #undef __FUNC__ 15775615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 15789b94acceSBarry Smith /*@ 15794b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 15809b94acceSBarry Smith 15819b94acceSBarry Smith Input Parameters: 15829b94acceSBarry Smith . snes - the SNES context 15839b94acceSBarry Smith . method - a known method 15849b94acceSBarry Smith 1585ae12b187SLois Curfman McInnes Options Database Command: 1586ae12b187SLois Curfman McInnes $ -snes_type <method> 1587ae12b187SLois Curfman McInnes $ Use -help for a list of available methods 1588ae12b187SLois Curfman McInnes $ (for instance, ls or tr) 1589ae12b187SLois Curfman McInnes 15909b94acceSBarry Smith Notes: 15919b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 15929b94acceSBarry Smith $ Systems of nonlinear equations: 159340191667SLois Curfman McInnes $ SNES_EQ_LS - Newton's method with line search 159440191667SLois Curfman McInnes $ SNES_EQ_TR - Newton's method with trust region 15959b94acceSBarry Smith $ Unconstrained minimization: 159640191667SLois Curfman McInnes $ SNES_UM_TR - Newton's method with trust region 159740191667SLois Curfman McInnes $ SNES_UM_LS - Newton's method with line search 15989b94acceSBarry Smith 1599ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1600ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1601ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1602ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1603ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1604ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1605ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1606ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1607ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1608ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1609a703fe33SLois Curfman McInnes 1610f536c699SSatish Balay .keywords: SNES, set, method 16119b94acceSBarry Smith @*/ 16124b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 16139b94acceSBarry Smith { 161484cb2905SBarry Smith int ierr; 16159b94acceSBarry Smith int (*r)(SNES); 16163a40ed3dSBarry Smith 16173a40ed3dSBarry Smith PetscFunctionBegin; 161877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16193a40ed3dSBarry Smith if (snes->type == (int) method) PetscFunctionReturn(0); 16207d1a2b2bSBarry Smith 16219b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 162284cb2905SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 162337fcc0dbSBarry Smith r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1624*a8c6a408SBarry Smith if (!r) {SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method");} 16250452661fSBarry Smith if (snes->data) PetscFree(snes->data); 16269b94acceSBarry Smith snes->set_method_called = 1; 16273a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 16283a40ed3dSBarry Smith PetscFunctionReturn(0); 16299b94acceSBarry Smith } 16309b94acceSBarry Smith 16319b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16325615d1e5SSatish Balay #undef __FUNC__ 1633d4bb536fSBarry Smith #define __FUNC__ "SNESRegister" 16349b94acceSBarry Smith /*@C 16359b94acceSBarry Smith SNESRegister - Adds the method to the nonlinear solver package, given 16364b0e389bSBarry Smith a function pointer and a nonlinear solver name of the type SNESType. 16379b94acceSBarry Smith 16389b94acceSBarry Smith Input Parameters: 16392d872ea7SLois Curfman McInnes . name - either a predefined name such as SNES_EQ_LS, or SNES_NEW 16402d872ea7SLois Curfman McInnes to indicate a new user-defined solver 164140191667SLois Curfman McInnes . sname - corresponding string for name 16429b94acceSBarry Smith . create - routine to create method context 16439b94acceSBarry Smith 164484cb2905SBarry Smith Output Parameter: 164584cb2905SBarry Smith . oname - type associated with this new solver 164684cb2905SBarry Smith 16472d872ea7SLois Curfman McInnes Notes: 16482d872ea7SLois Curfman McInnes Multiple user-defined nonlinear solvers can be added by calling 16492d872ea7SLois Curfman McInnes SNESRegister() with the input parameter "name" set to be SNES_NEW; 16502d872ea7SLois Curfman McInnes each call will return a unique solver type in the output 16512d872ea7SLois Curfman McInnes parameter "oname". 16522d872ea7SLois Curfman McInnes 16539b94acceSBarry Smith .keywords: SNES, nonlinear, register 16549b94acceSBarry Smith 16559b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy() 16569b94acceSBarry Smith @*/ 165784cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES)) 16589b94acceSBarry Smith { 16599b94acceSBarry Smith int ierr; 166084cb2905SBarry Smith static int numberregistered = 0; 166184cb2905SBarry Smith 16623a40ed3dSBarry Smith PetscFunctionBegin; 1663d252947aSBarry Smith if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++); 166484cb2905SBarry Smith 166584cb2905SBarry Smith if (oname) *oname = name; 166637fcc0dbSBarry Smith if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 166784cb2905SBarry Smith NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create ); 16683a40ed3dSBarry Smith PetscFunctionReturn(0); 16699b94acceSBarry Smith } 1670a847f771SSatish Balay 16719b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16725615d1e5SSatish Balay #undef __FUNC__ 1673d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy" 16749b94acceSBarry Smith /*@C 16759b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 16769b94acceSBarry Smith registered by SNESRegister(). 16779b94acceSBarry Smith 16789b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 16799b94acceSBarry Smith 16809b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 16819b94acceSBarry Smith @*/ 16829b94acceSBarry Smith int SNESRegisterDestroy() 16839b94acceSBarry Smith { 16843a40ed3dSBarry Smith PetscFunctionBegin; 168537fcc0dbSBarry Smith if (__SNESList) { 168637fcc0dbSBarry Smith NRDestroy( __SNESList ); 168737fcc0dbSBarry Smith __SNESList = 0; 16889b94acceSBarry Smith } 168984cb2905SBarry Smith SNESRegisterAllCalled = 0; 16903a40ed3dSBarry Smith PetscFunctionReturn(0); 16919b94acceSBarry Smith } 16929b94acceSBarry Smith 16935615d1e5SSatish Balay #undef __FUNC__ 1694d4bb536fSBarry Smith #define __FUNC__ "SNESGetType" 16959b94acceSBarry Smith /*@C 16969a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 16979b94acceSBarry Smith 16989b94acceSBarry Smith Input Parameter: 16994b0e389bSBarry Smith . snes - nonlinear solver context 17009b94acceSBarry Smith 17019b94acceSBarry Smith Output Parameter: 17029a28b0a6SLois Curfman McInnes . method - SNES method (or use PETSC_NULL) 17039a28b0a6SLois Curfman McInnes . name - name of SNES method (or use PETSC_NULL) 17049b94acceSBarry Smith 17059b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 17069b94acceSBarry Smith @*/ 17074b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name) 17089b94acceSBarry Smith { 170906a9b0adSLois Curfman McInnes int ierr; 17103a40ed3dSBarry Smith 17113a40ed3dSBarry Smith PetscFunctionBegin; 171237fcc0dbSBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 17134b0e389bSBarry Smith if (method) *method = (SNESType) snes->type; 171437fcc0dbSBarry Smith if (name) *name = NRFindName( __SNESList, (int) snes->type ); 17153a40ed3dSBarry Smith PetscFunctionReturn(0); 17169b94acceSBarry Smith } 17179b94acceSBarry Smith 17185615d1e5SSatish Balay #undef __FUNC__ 1719d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution" 17209b94acceSBarry Smith /*@C 17219b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 17229b94acceSBarry Smith stored. 17239b94acceSBarry Smith 17249b94acceSBarry Smith Input Parameter: 17259b94acceSBarry Smith . snes - the SNES context 17269b94acceSBarry Smith 17279b94acceSBarry Smith Output Parameter: 17289b94acceSBarry Smith . x - the solution 17299b94acceSBarry Smith 17309b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 17319b94acceSBarry Smith 1732a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 17339b94acceSBarry Smith @*/ 17349b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 17359b94acceSBarry Smith { 17363a40ed3dSBarry Smith PetscFunctionBegin; 173777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17389b94acceSBarry Smith *x = snes->vec_sol_always; 17393a40ed3dSBarry Smith PetscFunctionReturn(0); 17409b94acceSBarry Smith } 17419b94acceSBarry Smith 17425615d1e5SSatish Balay #undef __FUNC__ 1743d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" 17449b94acceSBarry Smith /*@C 17459b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 17469b94acceSBarry Smith stored. 17479b94acceSBarry Smith 17489b94acceSBarry Smith Input Parameter: 17499b94acceSBarry Smith . snes - the SNES context 17509b94acceSBarry Smith 17519b94acceSBarry Smith Output Parameter: 17529b94acceSBarry Smith . x - the solution update 17539b94acceSBarry Smith 17549b94acceSBarry Smith Notes: 17559b94acceSBarry Smith This vector is implementation dependent. 17569b94acceSBarry Smith 17579b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 17589b94acceSBarry Smith 17599b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 17609b94acceSBarry Smith @*/ 17619b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 17629b94acceSBarry Smith { 17633a40ed3dSBarry Smith PetscFunctionBegin; 176477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17659b94acceSBarry Smith *x = snes->vec_sol_update_always; 17663a40ed3dSBarry Smith PetscFunctionReturn(0); 17679b94acceSBarry Smith } 17689b94acceSBarry Smith 17695615d1e5SSatish Balay #undef __FUNC__ 1770d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction" 17719b94acceSBarry Smith /*@C 17723638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 17739b94acceSBarry Smith 17749b94acceSBarry Smith Input Parameter: 17759b94acceSBarry Smith . snes - the SNES context 17769b94acceSBarry Smith 17779b94acceSBarry Smith Output Parameter: 17783638b69dSLois Curfman McInnes . r - the function 17799b94acceSBarry Smith 17809b94acceSBarry Smith Notes: 17819b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 17829b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 17839b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 17849b94acceSBarry Smith 1785a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 17869b94acceSBarry Smith 178731615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 178831615d04SLois Curfman McInnes SNESGetGradient() 17899b94acceSBarry Smith @*/ 17909b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 17919b94acceSBarry Smith { 17923a40ed3dSBarry Smith PetscFunctionBegin; 179377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1794*a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1795*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1796*a8c6a408SBarry Smith } 17979b94acceSBarry Smith *r = snes->vec_func_always; 17983a40ed3dSBarry Smith PetscFunctionReturn(0); 17999b94acceSBarry Smith } 18009b94acceSBarry Smith 18015615d1e5SSatish Balay #undef __FUNC__ 1802d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient" 18039b94acceSBarry Smith /*@C 18043638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 18059b94acceSBarry Smith 18069b94acceSBarry Smith Input Parameter: 18079b94acceSBarry Smith . snes - the SNES context 18089b94acceSBarry Smith 18099b94acceSBarry Smith Output Parameter: 18109b94acceSBarry Smith . r - the gradient 18119b94acceSBarry Smith 18129b94acceSBarry Smith Notes: 18139b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 18149b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18159b94acceSBarry Smith SNESGetFunction(). 18169b94acceSBarry Smith 18179b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 18189b94acceSBarry Smith 181931615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 18209b94acceSBarry Smith @*/ 18219b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 18229b94acceSBarry Smith { 18233a40ed3dSBarry Smith PetscFunctionBegin; 182477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18253a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1826*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 18273a40ed3dSBarry Smith } 18289b94acceSBarry Smith *r = snes->vec_func_always; 18293a40ed3dSBarry Smith PetscFunctionReturn(0); 18309b94acceSBarry Smith } 18319b94acceSBarry Smith 18325615d1e5SSatish Balay #undef __FUNC__ 1833d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" 18349b94acceSBarry Smith /*@ 18359b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 18369b94acceSBarry Smith unconstrained minimization problems. 18379b94acceSBarry Smith 18389b94acceSBarry Smith Input Parameter: 18399b94acceSBarry Smith . snes - the SNES context 18409b94acceSBarry Smith 18419b94acceSBarry Smith Output Parameter: 18429b94acceSBarry Smith . r - the function 18439b94acceSBarry Smith 18449b94acceSBarry Smith Notes: 18459b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 18469b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18479b94acceSBarry Smith SNESGetFunction(). 18489b94acceSBarry Smith 18499b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 18509b94acceSBarry Smith 185131615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 18529b94acceSBarry Smith @*/ 18539b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 18549b94acceSBarry Smith { 18553a40ed3dSBarry Smith PetscFunctionBegin; 185677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 185774679c65SBarry Smith PetscValidScalarPointer(r); 18583a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1859*a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 18603a40ed3dSBarry Smith } 18619b94acceSBarry Smith *r = snes->fc; 18623a40ed3dSBarry Smith PetscFunctionReturn(0); 18639b94acceSBarry Smith } 18649b94acceSBarry Smith 18655615d1e5SSatish Balay #undef __FUNC__ 1866d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" 18673c7409f5SSatish Balay /*@C 18683c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1869d850072dSLois Curfman McInnes SNES options in the database. 18703c7409f5SSatish Balay 18713c7409f5SSatish Balay Input Parameter: 18723c7409f5SSatish Balay . snes - the SNES context 18733c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 18743c7409f5SSatish Balay 1875d850072dSLois Curfman McInnes Notes: 1876a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1877a83b1b31SSatish Balay The first character of all runtime options is AUTOMATICALLY the 1878a83b1b31SSatish Balay hyphen. 1879d850072dSLois Curfman McInnes 18803c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1881a86d99e1SLois Curfman McInnes 1882a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 18833c7409f5SSatish Balay @*/ 18843c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 18853c7409f5SSatish Balay { 18863c7409f5SSatish Balay int ierr; 18873c7409f5SSatish Balay 18883a40ed3dSBarry Smith PetscFunctionBegin; 188977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1890639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 18913c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 18923a40ed3dSBarry Smith PetscFunctionReturn(0); 18933c7409f5SSatish Balay } 18943c7409f5SSatish Balay 18955615d1e5SSatish Balay #undef __FUNC__ 1896d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" 18973c7409f5SSatish Balay /*@C 1898f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1899d850072dSLois Curfman McInnes SNES options in the database. 19003c7409f5SSatish Balay 19013c7409f5SSatish Balay Input Parameter: 19023c7409f5SSatish Balay . snes - the SNES context 19033c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 19043c7409f5SSatish Balay 1905d850072dSLois Curfman McInnes Notes: 1906a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1907a83b1b31SSatish Balay The first character of all runtime options is AUTOMATICALLY the 1908a83b1b31SSatish Balay hyphen. 1909d850072dSLois Curfman McInnes 19103c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1911a86d99e1SLois Curfman McInnes 1912a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 19133c7409f5SSatish Balay @*/ 19143c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 19153c7409f5SSatish Balay { 19163c7409f5SSatish Balay int ierr; 19173c7409f5SSatish Balay 19183a40ed3dSBarry Smith PetscFunctionBegin; 191977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1920639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19213c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19223a40ed3dSBarry Smith PetscFunctionReturn(0); 19233c7409f5SSatish Balay } 19243c7409f5SSatish Balay 19255615d1e5SSatish Balay #undef __FUNC__ 1926d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" 1927c83590e2SSatish Balay /*@ 19283c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 19293c7409f5SSatish Balay SNES options in the database. 19303c7409f5SSatish Balay 19313c7409f5SSatish Balay Input Parameter: 19323c7409f5SSatish Balay . snes - the SNES context 19333c7409f5SSatish Balay 19343c7409f5SSatish Balay Output Parameter: 19353c7409f5SSatish Balay . prefix - pointer to the prefix string used 19363c7409f5SSatish Balay 19373c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 1938a86d99e1SLois Curfman McInnes 1939a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 19403c7409f5SSatish Balay @*/ 19413c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 19423c7409f5SSatish Balay { 19433c7409f5SSatish Balay int ierr; 19443c7409f5SSatish Balay 19453a40ed3dSBarry Smith PetscFunctionBegin; 194677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1947639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19483a40ed3dSBarry Smith PetscFunctionReturn(0); 19493c7409f5SSatish Balay } 19503c7409f5SSatish Balay 1951ca161407SBarry Smith #undef __FUNC__ 1952ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp" 1953ca161407SBarry Smith /*@ 1954ca161407SBarry Smith SNESPrintHelp - Prints all options for the SNES component. 1955ca161407SBarry Smith 1956ca161407SBarry Smith Input Parameter: 1957ca161407SBarry Smith . snes - the SNES context 1958ca161407SBarry Smith 1959ca161407SBarry Smith Options Database Keys: 1960ca161407SBarry Smith $ -help, -h 1961ca161407SBarry Smith 1962ca161407SBarry Smith .keywords: SNES, nonlinear, help 1963ca161407SBarry Smith 1964ca161407SBarry Smith .seealso: SNESSetFromOptions() 1965ca161407SBarry Smith @*/ 1966ca161407SBarry Smith int SNESPrintHelp(SNES snes) 1967ca161407SBarry Smith { 1968ca161407SBarry Smith char p[64]; 1969ca161407SBarry Smith SNES_KSP_EW_ConvCtx *kctx; 1970ca161407SBarry Smith int ierr; 1971ca161407SBarry Smith 1972ca161407SBarry Smith PetscFunctionBegin; 1973ca161407SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1974ca161407SBarry Smith 1975ca161407SBarry Smith PetscStrcpy(p,"-"); 1976ca161407SBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 1977ca161407SBarry Smith 1978ca161407SBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 1979ca161407SBarry Smith 1980ca161407SBarry Smith PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 1981ca161407SBarry Smith ierr = NRPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",__SNESList);CHKERRQ(ierr); 1982ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 1983ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 1984ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 1985ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 1986ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 1987ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 1988ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 1989ca161407SBarry Smith PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 1990ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 1991ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 1992ca161407SBarry Smith residual norm at each iteration.\n",p); 1993ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 1994ca161407SBarry Smith residual norm for small residual norms. This is useful to conceal\n\ 1995ca161407SBarry Smith meaningless digits that may be different on different machines.\n",p); 1996ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 1997ca161407SBarry Smith if (snes->type == SNES_NONLINEAR_EQUATIONS) { 1998ca161407SBarry Smith PetscPrintf(snes->comm, 1999ca161407SBarry Smith " Options for solving systems of nonlinear equations only:\n"); 2000ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 2001ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 2002ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 2003ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p); 2004ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 2005ca161407SBarry Smith PetscPrintf(snes->comm, 2006ca161407SBarry Smith " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 2007ca161407SBarry Smith PetscPrintf(snes->comm, 2008ca161407SBarry Smith " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 2009ca161407SBarry Smith PetscPrintf(snes->comm, 2010ca161407SBarry Smith " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 2011ca161407SBarry Smith PetscPrintf(snes->comm, 2012ca161407SBarry Smith " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 2013ca161407SBarry Smith PetscPrintf(snes->comm, 2014ca161407SBarry Smith " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 2015ca161407SBarry Smith PetscPrintf(snes->comm, 2016ca161407SBarry Smith " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 2017ca161407SBarry Smith PetscPrintf(snes->comm, 2018ca161407SBarry Smith " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 2019ca161407SBarry Smith } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 2020ca161407SBarry Smith PetscPrintf(snes->comm, 2021ca161407SBarry Smith " Options for solving unconstrained minimization problems only:\n"); 2022ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 2023ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 2024ca161407SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 2025ca161407SBarry Smith } 2026ca161407SBarry Smith PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 2027ca161407SBarry Smith PetscPrintf(snes->comm,"a particular method\n"); 2028ca161407SBarry Smith if (snes->printhelp) { 2029ca161407SBarry Smith ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr); 2030ca161407SBarry Smith } 2031ca161407SBarry Smith PetscFunctionReturn(0); 2032ca161407SBarry Smith } 20333c7409f5SSatish Balay 20349b94acceSBarry Smith 20359b94acceSBarry Smith 20369b94acceSBarry Smith 20373a40ed3dSBarry Smith 2038