19b94acceSBarry Smith #ifndef lint 2*7a00f4a9SLois Curfman McInnes static char vcid[] = "$Id: snes.c,v 1.105 1997/01/06 20:29:45 balay Exp curfman $"; 39b94acceSBarry Smith #endif 49b94acceSBarry Smith 59b94acceSBarry Smith #include "draw.h" /*I "draw.h" I*/ 670f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 7f5eb4b81SSatish Balay #include "src/sys/nreg.h" 89b94acceSBarry Smith #include "pinclude/pviewer.h" 99b94acceSBarry Smith #include <math.h> 109b94acceSBarry Smith 11052efed2SBarry Smith extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*); 1233455573SSatish Balay extern int SNESPrintTypes_Private(MPI_Comm,char*,char*); 139b94acceSBarry Smith 145615d1e5SSatish Balay #undef __FUNC__ 155615d1e5SSatish Balay #define __FUNC__ "SNESView" 169b94acceSBarry Smith /*@ 179b94acceSBarry Smith SNESView - Prints the SNES data structure. 189b94acceSBarry Smith 199b94acceSBarry Smith Input Parameters: 209b94acceSBarry Smith . SNES - the SNES context 219b94acceSBarry Smith . viewer - visualization context 229b94acceSBarry Smith 239b94acceSBarry Smith Options Database Key: 249b94acceSBarry Smith $ -snes_view : calls SNESView() at end of SNESSolve() 259b94acceSBarry Smith 269b94acceSBarry Smith Notes: 279b94acceSBarry Smith The available visualization contexts include 286d4a8577SBarry Smith $ VIEWER_STDOUT_SELF - standard output (default) 296d4a8577SBarry Smith $ VIEWER_STDOUT_WORLD - synchronized standard 309b94acceSBarry Smith $ output where only the first processor opens 319b94acceSBarry Smith $ the file. All other processors send their 329b94acceSBarry Smith $ data to the first processor to print. 339b94acceSBarry Smith 349b94acceSBarry Smith The user can open alternative vistualization contexts with 35dbb450caSBarry Smith $ ViewerFileOpenASCII() - output to a specified file 369b94acceSBarry Smith 379b94acceSBarry Smith .keywords: SNES, view 389b94acceSBarry Smith 39dbb450caSBarry Smith .seealso: ViewerFileOpenASCII() 409b94acceSBarry Smith @*/ 419b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 429b94acceSBarry Smith { 439b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 449b94acceSBarry Smith FILE *fd; 459b94acceSBarry Smith int ierr; 469b94acceSBarry Smith SLES sles; 479b94acceSBarry Smith char *method; 4819bcc07fSBarry Smith ViewerType vtype; 499b94acceSBarry Smith 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); 67*7a00f4a9SLois Curfman McInnes PetscFPrintf(snes->comm,fd, 68*7a00f4a9SLois Curfman McInnes " total number of linear solver iterations=%d\n",snes->linear_its); 699b94acceSBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 7077c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 719b94acceSBarry Smith if (snes->ksp_ewconv) { 729b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 739b94acceSBarry Smith if (kctx) { 7477c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 759b94acceSBarry Smith " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 769b94acceSBarry Smith kctx->version); 7777c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 789b94acceSBarry Smith " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 799b94acceSBarry Smith kctx->rtol_max,kctx->threshold); 8077c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 819b94acceSBarry Smith kctx->gamma,kctx->alpha,kctx->alpha2); 829b94acceSBarry Smith } 839b94acceSBarry Smith } 84c30f285eSLois Curfman McInnes } else if (vtype == STRING_VIEWER) { 85c30f285eSLois Curfman McInnes SNESGetType(snes,PETSC_NULL,&method); 86c30f285eSLois Curfman McInnes ViewerStringSPrintf(viewer," %-3.3s",method); 8719bcc07fSBarry Smith } 889b94acceSBarry Smith SNESGetSLES(snes,&sles); 899b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 909b94acceSBarry Smith return 0; 919b94acceSBarry Smith } 929b94acceSBarry Smith 93639f9d9dSBarry Smith /* 94639f9d9dSBarry Smith We retain a list of functions that also take SNES command 95639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 96639f9d9dSBarry Smith */ 97639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 98639f9d9dSBarry Smith static int numberofsetfromoptions; 99639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 100639f9d9dSBarry Smith 1015615d1e5SSatish Balay #undef __FUNC__ 1025615d1e5SSatish Balay #define __FUNC__ "SNESAddOptionsChecker" 103639f9d9dSBarry Smith /*@ 104639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 105639f9d9dSBarry Smith 106639f9d9dSBarry Smith Input Parameter: 107639f9d9dSBarry Smith . snescheck - function that checks for options 108639f9d9dSBarry Smith 109639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 110639f9d9dSBarry Smith @*/ 111639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 112639f9d9dSBarry Smith { 113639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 114e3372554SBarry Smith SETERRQ(1,0,"Too many options checkers, only 5 allowed"); 115639f9d9dSBarry Smith } 116639f9d9dSBarry Smith 117639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 118639f9d9dSBarry Smith return 0; 119639f9d9dSBarry Smith } 120639f9d9dSBarry Smith 1215615d1e5SSatish Balay #undef __FUNC__ 1225615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1239b94acceSBarry Smith /*@ 124682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1259b94acceSBarry Smith 1269b94acceSBarry Smith Input Parameter: 1279b94acceSBarry Smith . snes - the SNES context 1289b94acceSBarry Smith 1299b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1309b94acceSBarry Smith 131a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1329b94acceSBarry Smith @*/ 1339b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1349b94acceSBarry Smith { 1354b0e389bSBarry Smith SNESType method; 1369b94acceSBarry Smith double tmp; 1379b94acceSBarry Smith SLES sles; 13881f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1393c7409f5SSatish Balay int version = PETSC_DEFAULT; 1409b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1419b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1429b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1439b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1449b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1459b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1469b94acceSBarry Smith 14781f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 14881f57ec7SBarry Smith 14981f57ec7SBarry Smith 15077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 151e3372554SBarry Smith if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp"); 152052efed2SBarry Smith ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 153052efed2SBarry Smith if (flg) { 154052efed2SBarry Smith ierr = SNESSetType(snes,method); CHKERRQ(ierr); 1559b94acceSBarry Smith } 1565334005bSBarry Smith else if (!snes->set_method_called) { 1575334005bSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 15840191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 1595334005bSBarry Smith } 1605334005bSBarry Smith else { 16140191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 1625334005bSBarry Smith } 1635334005bSBarry Smith } 1645334005bSBarry Smith 1653c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1663c7409f5SSatish Balay if (flg) { SNESPrintHelp(snes); } 1673c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 168d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);} 1693c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 170d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);} 1713c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 172d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); } 1733c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 1743c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 175d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 176d7a720efSLois Curfman McInnes if (flg) { SNESSetTrustRegionTolerance(snes,tmp); } 177d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 178d7a720efSLois Curfman McInnes if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); } 1793c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 1803c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 181b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 182b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 183b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 184b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 185b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 186b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 187b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 1889b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 1899b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 1905bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 1915bbad29bSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);} 1923c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 193639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 1943c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 195d31a3109SLois Curfman McInnes if (flg) { 196639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 197d31a3109SLois Curfman McInnes } 19881f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 1993c7409f5SSatish Balay if (flg) { 20017699dbbSLois Curfman McInnes int rank = 0; 201d7e8b826SBarry Smith DrawLG lg; 20217699dbbSLois Curfman McInnes MPI_Initialized(&rank); 20317699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 20417699dbbSLois Curfman McInnes if (!rank) { 20581f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2069b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 207f8c826e1SBarry Smith snes->xmonitor = lg; 2089b94acceSBarry Smith PLogObjectParent(snes,lg); 2099b94acceSBarry Smith } 2109b94acceSBarry Smith } 2113c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2123c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2139b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2149b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 21594a424c1SBarry Smith PLogInfo(snes, 21631615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 21731615d04SLois Curfman McInnes } 21831615d04SLois Curfman McInnes else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 21931615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 22031615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 22194a424c1SBarry Smith PLogInfo(snes, 22231615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2239b94acceSBarry Smith } 224639f9d9dSBarry Smith 225639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 226639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 227639f9d9dSBarry Smith } 228639f9d9dSBarry Smith 2299b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2309b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 2319b94acceSBarry Smith if (!snes->setfromoptions) return 0; 2329b94acceSBarry Smith return (*snes->setfromoptions)(snes); 2339b94acceSBarry Smith } 2349b94acceSBarry Smith 2355615d1e5SSatish Balay #undef __FUNC__ 2365615d1e5SSatish Balay #define __FUNC__ "SNESPrintHelp" 2379b94acceSBarry Smith /*@ 2389b94acceSBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2399b94acceSBarry Smith 2409b94acceSBarry Smith Input Parameter: 2419b94acceSBarry Smith . snes - the SNES context 2429b94acceSBarry Smith 243a703fe33SLois Curfman McInnes Options Database Keys: 244a703fe33SLois Curfman McInnes $ -help, -h 245a703fe33SLois Curfman McInnes 2469b94acceSBarry Smith .keywords: SNES, nonlinear, help 2479b94acceSBarry Smith 248682d7d0cSBarry Smith .seealso: SNESSetFromOptions() 2499b94acceSBarry Smith @*/ 2509b94acceSBarry Smith int SNESPrintHelp(SNES snes) 2519b94acceSBarry Smith { 2526daaf66cSBarry Smith char p[64]; 2536daaf66cSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2546daaf66cSBarry Smith 25577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2566daaf66cSBarry Smith 2576daaf66cSBarry Smith PetscStrcpy(p,"-"); 2586daaf66cSBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2596daaf66cSBarry Smith 2606daaf66cSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2616daaf66cSBarry Smith 262d31a3109SLois Curfman McInnes PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 26333455573SSatish Balay SNESPrintTypes_Private(snes->comm,p,"snes_type"); 26477c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 265b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 266b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 267b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 268b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 269b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 270b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 2715bbad29bSBarry Smith PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 2725bbad29bSBarry Smith PetscPrintf(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 273d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor\n",p); 274d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 275b18e04deSLois Curfman McInnes if (snes->type == SNES_NONLINEAR_EQUATIONS) { 27677c4ece6SBarry Smith PetscPrintf(snes->comm, 277d31a3109SLois Curfman McInnes " Options for solving systems of nonlinear equations only:\n"); 27877c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 27977c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 280ef1dfb11SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 2811650c7b0SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p); 2821650c7b0SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p); 28377c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 28477c4ece6SBarry Smith PetscPrintf(snes->comm, 285b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 28677c4ece6SBarry Smith PetscPrintf(snes->comm, 287b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 28877c4ece6SBarry Smith PetscPrintf(snes->comm, 289b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 29077c4ece6SBarry Smith PetscPrintf(snes->comm, 291b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 29277c4ece6SBarry Smith PetscPrintf(snes->comm, 293b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 29477c4ece6SBarry Smith PetscPrintf(snes->comm, 295b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 29677c4ece6SBarry Smith PetscPrintf(snes->comm, 297b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 298b18e04deSLois Curfman McInnes } 299b18e04deSLois Curfman McInnes else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 30077c4ece6SBarry Smith PetscPrintf(snes->comm, 301d31a3109SLois Curfman McInnes " Options for solving unconstrained minimization problems only:\n"); 302b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 30377c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 30477c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 305d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p); 306d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p); 307b18e04deSLois Curfman McInnes } 3084537a946SLois Curfman McInnes PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 30977c4ece6SBarry Smith PetscPrintf(snes->comm,"a particular method\n"); 3106daaf66cSBarry Smith if (snes->printhelp) (*snes->printhelp)(snes,p); 3119b94acceSBarry Smith return 0; 3129b94acceSBarry Smith } 313a847f771SSatish Balay 3145615d1e5SSatish Balay #undef __FUNC__ 3155615d1e5SSatish Balay #define __FUNC__ "SNESSetApplicationContext" 3169b94acceSBarry Smith /*@ 3179b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 3189b94acceSBarry Smith the nonlinear solvers. 3199b94acceSBarry Smith 3209b94acceSBarry Smith Input Parameters: 3219b94acceSBarry Smith . snes - the SNES context 3229b94acceSBarry Smith . usrP - optional user context 3239b94acceSBarry Smith 3249b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3259b94acceSBarry Smith 3269b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3279b94acceSBarry Smith @*/ 3289b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3299b94acceSBarry Smith { 33077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3319b94acceSBarry Smith snes->user = usrP; 3329b94acceSBarry Smith return 0; 3339b94acceSBarry Smith } 33474679c65SBarry Smith 3355615d1e5SSatish Balay #undef __FUNC__ 3365615d1e5SSatish Balay #define __FUNC__ "SNESGetApplicationContext" 3379b94acceSBarry Smith /*@C 3389b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3399b94acceSBarry Smith nonlinear solvers. 3409b94acceSBarry Smith 3419b94acceSBarry Smith Input Parameter: 3429b94acceSBarry Smith . snes - SNES context 3439b94acceSBarry Smith 3449b94acceSBarry Smith Output Parameter: 3459b94acceSBarry Smith . usrP - user context 3469b94acceSBarry Smith 3479b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3489b94acceSBarry Smith 3499b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3509b94acceSBarry Smith @*/ 3519b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3529b94acceSBarry Smith { 35377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3549b94acceSBarry Smith *usrP = snes->user; 3559b94acceSBarry Smith return 0; 3569b94acceSBarry Smith } 35774679c65SBarry Smith 3585615d1e5SSatish Balay #undef __FUNC__ 3595615d1e5SSatish Balay #define __FUNC__ "SNESGetIterationNumber" 3609b94acceSBarry Smith /*@ 3619b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3629b94acceSBarry Smith nonlinear solver. 3639b94acceSBarry Smith 3649b94acceSBarry Smith Input Parameter: 3659b94acceSBarry Smith . snes - SNES context 3669b94acceSBarry Smith 3679b94acceSBarry Smith Output Parameter: 3689b94acceSBarry Smith . iter - iteration number 3699b94acceSBarry Smith 3709b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3719b94acceSBarry Smith @*/ 3729b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3739b94acceSBarry Smith { 37477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 37574679c65SBarry Smith PetscValidIntPointer(iter); 3769b94acceSBarry Smith *iter = snes->iter; 3779b94acceSBarry Smith return 0; 3789b94acceSBarry Smith } 37974679c65SBarry Smith 3805615d1e5SSatish Balay #undef __FUNC__ 3815615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 3829b94acceSBarry Smith /*@ 3839b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3849b94acceSBarry Smith with SNESSSetFunction(). 3859b94acceSBarry Smith 3869b94acceSBarry Smith Input Parameter: 3879b94acceSBarry Smith . snes - SNES context 3889b94acceSBarry Smith 3899b94acceSBarry Smith Output Parameter: 3909b94acceSBarry Smith . fnorm - 2-norm of function 3919b94acceSBarry Smith 3929b94acceSBarry Smith Note: 3939b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 394a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 395a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 3969b94acceSBarry Smith 3979b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 398a86d99e1SLois Curfman McInnes 399a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 4009b94acceSBarry Smith @*/ 4019b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 4029b94acceSBarry Smith { 40377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 40474679c65SBarry Smith PetscValidScalarPointer(fnorm); 40574679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 406e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 40774679c65SBarry Smith } 4089b94acceSBarry Smith *fnorm = snes->norm; 4099b94acceSBarry Smith return 0; 4109b94acceSBarry Smith } 41174679c65SBarry Smith 4125615d1e5SSatish Balay #undef __FUNC__ 4135615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 4149b94acceSBarry Smith /*@ 4159b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 4169b94acceSBarry Smith with SNESSSetGradient(). 4179b94acceSBarry Smith 4189b94acceSBarry Smith Input Parameter: 4199b94acceSBarry Smith . snes - SNES context 4209b94acceSBarry Smith 4219b94acceSBarry Smith Output Parameter: 4229b94acceSBarry Smith . fnorm - 2-norm of gradient 4239b94acceSBarry Smith 4249b94acceSBarry Smith Note: 4259b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 426a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 427a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4289b94acceSBarry Smith 4299b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 430a86d99e1SLois Curfman McInnes 431a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4329b94acceSBarry Smith @*/ 4339b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4349b94acceSBarry Smith { 43577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 43674679c65SBarry Smith PetscValidScalarPointer(gnorm); 43774679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 438e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 43974679c65SBarry Smith } 4409b94acceSBarry Smith *gnorm = snes->norm; 4419b94acceSBarry Smith return 0; 4429b94acceSBarry Smith } 44374679c65SBarry Smith 4445615d1e5SSatish Balay #undef __FUNC__ 4455615d1e5SSatish Balay #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 4469b94acceSBarry Smith /*@ 4479b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4489b94acceSBarry Smith attempted by the nonlinear solver. 4499b94acceSBarry Smith 4509b94acceSBarry Smith Input Parameter: 4519b94acceSBarry Smith . snes - SNES context 4529b94acceSBarry Smith 4539b94acceSBarry Smith Output Parameter: 4549b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4559b94acceSBarry Smith 4569b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4579b94acceSBarry Smith @*/ 4589b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4599b94acceSBarry Smith { 46077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 46174679c65SBarry Smith PetscValidIntPointer(nfails); 4629b94acceSBarry Smith *nfails = snes->nfailures; 4639b94acceSBarry Smith return 0; 4649b94acceSBarry Smith } 465a847f771SSatish Balay 4665615d1e5SSatish Balay #undef __FUNC__ 4675615d1e5SSatish Balay #define __FUNC__ "SNESGetSLES" 4689b94acceSBarry Smith /*@C 4699b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 4709b94acceSBarry Smith 4719b94acceSBarry Smith Input Parameter: 4729b94acceSBarry Smith . snes - the SNES context 4739b94acceSBarry Smith 4749b94acceSBarry Smith Output Parameter: 4759b94acceSBarry Smith . sles - the SLES context 4769b94acceSBarry Smith 4779b94acceSBarry Smith Notes: 4789b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 4799b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 4809b94acceSBarry Smith KSP and PC contexts as well. 4819b94acceSBarry Smith 4829b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 4839b94acceSBarry Smith 4849b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 4859b94acceSBarry Smith @*/ 4869b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 4879b94acceSBarry Smith { 48877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 4899b94acceSBarry Smith *sles = snes->sles; 4909b94acceSBarry Smith return 0; 4919b94acceSBarry Smith } 4929b94acceSBarry Smith /* -----------------------------------------------------------*/ 4935615d1e5SSatish Balay #undef __FUNC__ 4945615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 4959b94acceSBarry Smith /*@C 4969b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 4979b94acceSBarry Smith 4989b94acceSBarry Smith Input Parameter: 4999b94acceSBarry Smith . comm - MPI communicator 5009b94acceSBarry Smith . type - type of method, one of 5019b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS 5029b94acceSBarry Smith $ (for systems of nonlinear equations) 5039b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION 5049b94acceSBarry Smith $ (for unconstrained minimization) 5059b94acceSBarry Smith 5069b94acceSBarry Smith Output Parameter: 5079b94acceSBarry Smith . outsnes - the new SNES context 5089b94acceSBarry Smith 509c1f60f51SBarry Smith Options Database Key: 51019bd6747SLois Curfman McInnes $ -snes_mf - use default matrix-free Jacobian-vector products, 51119bd6747SLois Curfman McInnes $ and no preconditioning matrix 51219bd6747SLois Curfman McInnes $ -snes_mf_operator - use default matrix-free Jacobian-vector 51319bd6747SLois Curfman McInnes $ products, and a user-provided preconditioning matrix 51419bd6747SLois Curfman McInnes $ as set by SNESSetJacobian() 51519bd6747SLois Curfman McInnes $ -snes_fd - use (slow!) finite differences to compute Jacobian 516c1f60f51SBarry Smith 5179b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5189b94acceSBarry Smith 51963a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 5209b94acceSBarry Smith @*/ 5214b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 5229b94acceSBarry Smith { 5239b94acceSBarry Smith int ierr; 5249b94acceSBarry Smith SNES snes; 5259b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 52637fcc0dbSBarry Smith 5279b94acceSBarry Smith *outsnes = 0; 5282a0acf01SLois Curfman McInnes if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 529e3372554SBarry Smith SETERRQ(1,0,"incorrect method type"); 5300452661fSBarry Smith PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm); 5319b94acceSBarry Smith PLogObjectCreate(snes); 5329b94acceSBarry Smith snes->max_its = 50; 5339b94acceSBarry Smith snes->max_funcs = 1000; 5349b94acceSBarry Smith snes->norm = 0.0; 5355a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 536b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 537b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5389b94acceSBarry Smith snes->atol = 1.e-10; 5395a2d0531SBarry Smith } 540b4874afaSBarry Smith else { 541b4874afaSBarry Smith snes->rtol = 1.e-8; 542b4874afaSBarry Smith snes->ttol = 0.0; 543b4874afaSBarry Smith snes->atol = 1.e-50; 544b4874afaSBarry Smith } 5459b94acceSBarry Smith snes->xtol = 1.e-8; 546b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5479b94acceSBarry Smith snes->nfuncs = 0; 5489b94acceSBarry Smith snes->nfailures = 0; 549*7a00f4a9SLois Curfman McInnes snes->linear_its = 0; 550639f9d9dSBarry Smith snes->numbermonitors = 0; 5519b94acceSBarry Smith snes->data = 0; 5529b94acceSBarry Smith snes->view = 0; 5539b94acceSBarry Smith snes->computeumfunction = 0; 5549b94acceSBarry Smith snes->umfunP = 0; 5559b94acceSBarry Smith snes->fc = 0; 5569b94acceSBarry Smith snes->deltatol = 1.e-12; 5579b94acceSBarry Smith snes->fmin = -1.e30; 5589b94acceSBarry Smith snes->method_class = type; 5599b94acceSBarry Smith snes->set_method_called = 0; 560a703fe33SLois Curfman McInnes snes->setup_called = 0; 5619b94acceSBarry Smith snes->ksp_ewconv = 0; 5627d1a2b2bSBarry Smith snes->type = -1; 5636f24a144SLois Curfman McInnes snes->xmonitor = 0; 5646f24a144SLois Curfman McInnes snes->vwork = 0; 5656f24a144SLois Curfman McInnes snes->nwork = 0; 5669b94acceSBarry Smith 5679b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 5680452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 5699b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 5709b94acceSBarry Smith kctx->version = 2; 5719b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 5729b94acceSBarry Smith this was too large for some test cases */ 5739b94acceSBarry Smith kctx->rtol_last = 0; 5749b94acceSBarry Smith kctx->rtol_max = .9; 5759b94acceSBarry Smith kctx->gamma = 1.0; 5769b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 5779b94acceSBarry Smith kctx->alpha = kctx->alpha2; 5789b94acceSBarry Smith kctx->threshold = .1; 5799b94acceSBarry Smith kctx->lresid_last = 0; 5809b94acceSBarry Smith kctx->norm_last = 0; 5819b94acceSBarry Smith 5829b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 5839b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 5845334005bSBarry Smith 5859b94acceSBarry Smith *outsnes = snes; 5869b94acceSBarry Smith return 0; 5879b94acceSBarry Smith } 5889b94acceSBarry Smith 5899b94acceSBarry Smith /* --------------------------------------------------------------- */ 5905615d1e5SSatish Balay #undef __FUNC__ 5915615d1e5SSatish Balay #define __FUNC__ "SNESSetFunction" 5929b94acceSBarry Smith /*@C 5939b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 5949b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 5959b94acceSBarry Smith equations. 5969b94acceSBarry Smith 5979b94acceSBarry Smith Input Parameters: 5989b94acceSBarry Smith . snes - the SNES context 5999b94acceSBarry Smith . func - function evaluation routine 6009b94acceSBarry Smith . r - vector to store function value 6012cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 6022cd2dfdcSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6039b94acceSBarry Smith 6049b94acceSBarry Smith Calling sequence of func: 605313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Vec f,void *ctx); 6069b94acceSBarry Smith 6079b94acceSBarry Smith . x - input vector 608313e4042SLois Curfman McInnes . f - function vector 6092cd2dfdcSLois Curfman McInnes . ctx - optional user-defined function context 6109b94acceSBarry Smith 6119b94acceSBarry Smith Notes: 6129b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6139b94acceSBarry Smith $ f'(x) x = -f(x), 6149b94acceSBarry Smith $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 6159b94acceSBarry Smith 6169b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6179b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6189b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 6199b94acceSBarry Smith 6209b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6219b94acceSBarry Smith 622a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6239b94acceSBarry Smith @*/ 6245334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 6259b94acceSBarry Smith { 62677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 627e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6289b94acceSBarry Smith snes->computefunction = func; 6299b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6309b94acceSBarry Smith snes->funP = ctx; 6319b94acceSBarry Smith return 0; 6329b94acceSBarry Smith } 6339b94acceSBarry Smith 6345615d1e5SSatish Balay #undef __FUNC__ 6355615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 6369b94acceSBarry Smith /*@ 6379b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6389b94acceSBarry Smith SNESSetFunction(). 6399b94acceSBarry Smith 6409b94acceSBarry Smith Input Parameters: 6419b94acceSBarry Smith . snes - the SNES context 6429b94acceSBarry Smith . x - input vector 6439b94acceSBarry Smith 6449b94acceSBarry Smith Output Parameter: 6453638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6469b94acceSBarry Smith 6471bffabb2SLois Curfman McInnes Notes: 6489b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6499b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6509b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6519b94acceSBarry Smith 6529b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6539b94acceSBarry Smith 654a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6559b94acceSBarry Smith @*/ 6569b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 6579b94acceSBarry Smith { 6589b94acceSBarry Smith int ierr; 6599b94acceSBarry Smith 66074679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 661e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6629b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 6639b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 6649b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 6659b94acceSBarry Smith return 0; 6669b94acceSBarry Smith } 6679b94acceSBarry Smith 6685615d1e5SSatish Balay #undef __FUNC__ 6695615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunction" 6709b94acceSBarry Smith /*@C 6719b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 6729b94acceSBarry Smith unconstrained minimization. 6739b94acceSBarry Smith 6749b94acceSBarry Smith Input Parameters: 6759b94acceSBarry Smith . snes - the SNES context 6769b94acceSBarry Smith . func - function evaluation routine 677044dda88SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 678044dda88SLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6799b94acceSBarry Smith 6809b94acceSBarry Smith Calling sequence of func: 6819b94acceSBarry Smith . func (SNES snes,Vec x,double *f,void *ctx); 6829b94acceSBarry Smith 6839b94acceSBarry Smith . x - input vector 6849b94acceSBarry Smith . f - function 685044dda88SLois Curfman McInnes . ctx - [optional] user-defined function context 6869b94acceSBarry Smith 6879b94acceSBarry Smith Notes: 6889b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 6899b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 6909b94acceSBarry Smith SNESSetFunction(). 6919b94acceSBarry Smith 6929b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 6939b94acceSBarry Smith 694a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 695a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 6969b94acceSBarry Smith @*/ 6979b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 6989b94acceSBarry Smith void *ctx) 6999b94acceSBarry Smith { 70077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 701e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7029b94acceSBarry Smith snes->computeumfunction = func; 7039b94acceSBarry Smith snes->umfunP = ctx; 7049b94acceSBarry Smith return 0; 7059b94acceSBarry Smith } 7069b94acceSBarry Smith 7075615d1e5SSatish Balay #undef __FUNC__ 7085615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 7099b94acceSBarry Smith /*@ 7109b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 7119b94acceSBarry Smith set with SNESSetMinimizationFunction(). 7129b94acceSBarry Smith 7139b94acceSBarry Smith Input Parameters: 7149b94acceSBarry Smith . snes - the SNES context 7159b94acceSBarry Smith . x - input vector 7169b94acceSBarry Smith 7179b94acceSBarry Smith Output Parameter: 7189b94acceSBarry Smith . y - function value 7199b94acceSBarry Smith 7209b94acceSBarry Smith Notes: 7219b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 7229b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7239b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 724a86d99e1SLois Curfman McInnes 725a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 726a86d99e1SLois Curfman McInnes 727a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 728a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 7299b94acceSBarry Smith @*/ 7309b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 7319b94acceSBarry Smith { 7329b94acceSBarry Smith int ierr; 733e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7349b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 7359b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 7369b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7379b94acceSBarry Smith return 0; 7389b94acceSBarry Smith } 7399b94acceSBarry Smith 7405615d1e5SSatish Balay #undef __FUNC__ 7415615d1e5SSatish Balay #define __FUNC__ "SNESSetGradient" 7429b94acceSBarry Smith /*@C 7439b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 7449b94acceSBarry Smith vector for use by the SNES routines. 7459b94acceSBarry Smith 7469b94acceSBarry Smith Input Parameters: 7479b94acceSBarry Smith . snes - the SNES context 7489b94acceSBarry Smith . func - function evaluation routine 749044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 750044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 7519b94acceSBarry Smith . r - vector to store gradient value 7529b94acceSBarry Smith 7539b94acceSBarry Smith Calling sequence of func: 7549b94acceSBarry Smith . func (SNES, Vec x, Vec g, void *ctx); 7559b94acceSBarry Smith 7569b94acceSBarry Smith . x - input vector 7579b94acceSBarry Smith . g - gradient vector 758044dda88SLois Curfman McInnes . ctx - optional user-defined gradient context 7599b94acceSBarry Smith 7609b94acceSBarry Smith Notes: 7619b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7629b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7639b94acceSBarry Smith SNESSetFunction(). 7649b94acceSBarry Smith 7659b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 7669b94acceSBarry Smith 767a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 768a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 7699b94acceSBarry Smith @*/ 77074679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 7719b94acceSBarry Smith { 77277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 773e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 7749b94acceSBarry Smith snes->computefunction = func; 7759b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7769b94acceSBarry Smith snes->funP = ctx; 7779b94acceSBarry Smith return 0; 7789b94acceSBarry Smith } 7799b94acceSBarry Smith 7805615d1e5SSatish Balay #undef __FUNC__ 7815615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 7829b94acceSBarry Smith /*@ 783a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 784a86d99e1SLois Curfman McInnes SNESSetGradient(). 7859b94acceSBarry Smith 7869b94acceSBarry Smith Input Parameters: 7879b94acceSBarry Smith . snes - the SNES context 7889b94acceSBarry Smith . x - input vector 7899b94acceSBarry Smith 7909b94acceSBarry Smith Output Parameter: 7919b94acceSBarry Smith . y - gradient vector 7929b94acceSBarry Smith 7939b94acceSBarry Smith Notes: 7949b94acceSBarry Smith SNESComputeGradient() is valid only for 7959b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7969b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 797a86d99e1SLois Curfman McInnes 798a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 799a86d99e1SLois Curfman McInnes 800a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 801a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 8029b94acceSBarry Smith @*/ 8039b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 8049b94acceSBarry Smith { 8059b94acceSBarry Smith int ierr; 80674679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 807e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8089b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 8099b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 8109b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 8119b94acceSBarry Smith return 0; 8129b94acceSBarry Smith } 8139b94acceSBarry Smith 8145615d1e5SSatish Balay #undef __FUNC__ 8155615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 81662fef451SLois Curfman McInnes /*@ 81762fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 81862fef451SLois Curfman McInnes set with SNESSetJacobian(). 81962fef451SLois Curfman McInnes 82062fef451SLois Curfman McInnes Input Parameters: 82162fef451SLois Curfman McInnes . snes - the SNES context 82262fef451SLois Curfman McInnes . x - input vector 82362fef451SLois Curfman McInnes 82462fef451SLois Curfman McInnes Output Parameters: 82562fef451SLois Curfman McInnes . A - Jacobian matrix 82662fef451SLois Curfman McInnes . B - optional preconditioning matrix 82762fef451SLois Curfman McInnes . flag - flag indicating matrix structure 82862fef451SLois Curfman McInnes 82962fef451SLois Curfman McInnes Notes: 83062fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 83162fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 83262fef451SLois Curfman McInnes 833dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 834dc5a77f8SLois Curfman McInnes flag parameter. 83562fef451SLois Curfman McInnes 83662fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 83762fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 83862fef451SLois Curfman McInnes methods is SNESComputeJacobian(). 83962fef451SLois Curfman McInnes 84062fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 84162fef451SLois Curfman McInnes 84262fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 84362fef451SLois Curfman McInnes @*/ 8449b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8459b94acceSBarry Smith { 8469b94acceSBarry Smith int ierr; 84774679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 848e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 8499b94acceSBarry Smith if (!snes->computejacobian) return 0; 8509b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 851c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 8529b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 8539b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 8546d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 85577c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 85677c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 8579b94acceSBarry Smith return 0; 8589b94acceSBarry Smith } 8599b94acceSBarry Smith 8605615d1e5SSatish Balay #undef __FUNC__ 8615615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 86262fef451SLois Curfman McInnes /*@ 86362fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 86462fef451SLois Curfman McInnes set with SNESSetHessian(). 86562fef451SLois Curfman McInnes 86662fef451SLois Curfman McInnes Input Parameters: 86762fef451SLois Curfman McInnes . snes - the SNES context 86862fef451SLois Curfman McInnes . x - input vector 86962fef451SLois Curfman McInnes 87062fef451SLois Curfman McInnes Output Parameters: 87162fef451SLois Curfman McInnes . A - Hessian matrix 87262fef451SLois Curfman McInnes . B - optional preconditioning matrix 87362fef451SLois Curfman McInnes . flag - flag indicating matrix structure 87462fef451SLois Curfman McInnes 87562fef451SLois Curfman McInnes Notes: 87662fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 87762fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 87862fef451SLois Curfman McInnes 879dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 880dc5a77f8SLois Curfman McInnes flag parameter. 88162fef451SLois Curfman McInnes 88262fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 88362fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 88462fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 88562fef451SLois Curfman McInnes 88662fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 88762fef451SLois Curfman McInnes 888a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 889a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 89062fef451SLois Curfman McInnes @*/ 89162fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 8929b94acceSBarry Smith { 8939b94acceSBarry Smith int ierr; 89474679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 895e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8969b94acceSBarry Smith if (!snes->computejacobian) return 0; 89762fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 8980f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 89962fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 90062fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 9010f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 90277c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 90377c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9049b94acceSBarry Smith return 0; 9059b94acceSBarry Smith } 9069b94acceSBarry Smith 9075615d1e5SSatish Balay #undef __FUNC__ 9085615d1e5SSatish Balay #define __FUNC__ "SNESSetJacobian" 9099b94acceSBarry Smith /*@C 9109b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 911044dda88SLois Curfman McInnes location to store the matrix. 9129b94acceSBarry Smith 9139b94acceSBarry Smith Input Parameters: 9149b94acceSBarry Smith . snes - the SNES context 9159b94acceSBarry Smith . A - Jacobian matrix 9169b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 9179b94acceSBarry Smith . func - Jacobian evaluation routine 9182cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 9192cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 9209b94acceSBarry Smith 9219b94acceSBarry Smith Calling sequence of func: 922313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9239b94acceSBarry Smith 9249b94acceSBarry Smith . x - input vector 9259b94acceSBarry Smith . A - Jacobian matrix 9269b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 927ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 928ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 9292cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Jacobian context 9309b94acceSBarry Smith 9319b94acceSBarry Smith Notes: 932dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 9332cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 934ac21db08SLois Curfman McInnes 935ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9369b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9379b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9389b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9399b94acceSBarry Smith throughout the global iterations. 9409b94acceSBarry Smith 9419b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9429b94acceSBarry Smith 943ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 9449b94acceSBarry Smith @*/ 9459b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 9469b94acceSBarry Smith MatStructure*,void*),void *ctx) 9479b94acceSBarry Smith { 94877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 94974679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 950e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 9519b94acceSBarry Smith snes->computejacobian = func; 9529b94acceSBarry Smith snes->jacP = ctx; 9539b94acceSBarry Smith snes->jacobian = A; 9549b94acceSBarry Smith snes->jacobian_pre = B; 9559b94acceSBarry Smith return 0; 9569b94acceSBarry Smith } 95762fef451SLois Curfman McInnes 9585615d1e5SSatish Balay #undef __FUNC__ 9595615d1e5SSatish Balay #define __FUNC__ "SNESGetJacobian" 960b4fd4287SBarry Smith /*@ 961b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 962b4fd4287SBarry Smith provided context for evaluating the Jacobian. 963b4fd4287SBarry Smith 964b4fd4287SBarry Smith Input Parameter: 965b4fd4287SBarry Smith . snes - the nonlinear solver context 966b4fd4287SBarry Smith 967b4fd4287SBarry Smith Output Parameters: 968b4fd4287SBarry Smith . A - location to stash Jacobian matrix (or PETSC_NULL) 969b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 970b4fd4287SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 971b4fd4287SBarry Smith 972b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 973b4fd4287SBarry Smith @*/ 974b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 975b4fd4287SBarry Smith { 97674679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 977e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 978b4fd4287SBarry Smith if (A) *A = snes->jacobian; 979b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 980b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 981b4fd4287SBarry Smith return 0; 982b4fd4287SBarry Smith } 983b4fd4287SBarry Smith 9845615d1e5SSatish Balay #undef __FUNC__ 9855615d1e5SSatish Balay #define __FUNC__ "SNESSetHessian" 9869b94acceSBarry Smith /*@C 9879b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 988044dda88SLois Curfman McInnes location to store the matrix. 9899b94acceSBarry Smith 9909b94acceSBarry Smith Input Parameters: 9919b94acceSBarry Smith . snes - the SNES context 9929b94acceSBarry Smith . A - Hessian matrix 9939b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 9949b94acceSBarry Smith . func - Jacobian evaluation routine 995313e4042SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 996313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 9979b94acceSBarry Smith 9989b94acceSBarry Smith Calling sequence of func: 999313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 10009b94acceSBarry Smith 10019b94acceSBarry Smith . x - input vector 10029b94acceSBarry Smith . A - Hessian matrix 10039b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1004ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1005ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 10062cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Hessian context 10079b94acceSBarry Smith 10089b94acceSBarry Smith Notes: 1009dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10102cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1011ac21db08SLois Curfman McInnes 10129b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 10139b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 10149b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10159b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10169b94acceSBarry Smith throughout the global iterations. 10179b94acceSBarry Smith 10189b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 10199b94acceSBarry Smith 1020ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 10219b94acceSBarry Smith @*/ 10229b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 10239b94acceSBarry Smith MatStructure*,void*),void *ctx) 10249b94acceSBarry Smith { 102577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 102674679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1027e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 10289b94acceSBarry Smith snes->computejacobian = func; 10299b94acceSBarry Smith snes->jacP = ctx; 10309b94acceSBarry Smith snes->jacobian = A; 10319b94acceSBarry Smith snes->jacobian_pre = B; 10329b94acceSBarry Smith return 0; 10339b94acceSBarry Smith } 10349b94acceSBarry Smith 10355615d1e5SSatish Balay #undef __FUNC__ 10365615d1e5SSatish Balay #define __FUNC__ "SNESGetHessian" 103762fef451SLois Curfman McInnes /*@ 103862fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 103962fef451SLois Curfman McInnes provided context for evaluating the Hessian. 104062fef451SLois Curfman McInnes 104162fef451SLois Curfman McInnes Input Parameter: 104262fef451SLois Curfman McInnes . snes - the nonlinear solver context 104362fef451SLois Curfman McInnes 104462fef451SLois Curfman McInnes Output Parameters: 104562fef451SLois Curfman McInnes . A - location to stash Hessian matrix (or PETSC_NULL) 104662fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 104762fef451SLois Curfman McInnes . ctx - location to stash Hessian ctx (or PETSC_NULL) 104862fef451SLois Curfman McInnes 104962fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 105062fef451SLois Curfman McInnes @*/ 105162fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 105262fef451SLois Curfman McInnes { 105374679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1054e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 105562fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 105662fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 105762fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 105862fef451SLois Curfman McInnes return 0; 105962fef451SLois Curfman McInnes } 106062fef451SLois Curfman McInnes 10619b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 10629b94acceSBarry Smith 10635615d1e5SSatish Balay #undef __FUNC__ 10645615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 10659b94acceSBarry Smith /*@ 10669b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1067272ac6f2SLois Curfman McInnes of a nonlinear solver. 10689b94acceSBarry Smith 10699b94acceSBarry Smith Input Parameter: 10709b94acceSBarry Smith . snes - the SNES context 10718ddd3da0SLois Curfman McInnes . x - the solution vector 10729b94acceSBarry Smith 1073272ac6f2SLois Curfman McInnes Notes: 1074272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1075272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1076272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1077272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1078272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1079272ac6f2SLois Curfman McInnes 10809b94acceSBarry Smith .keywords: SNES, nonlinear, setup 10819b94acceSBarry Smith 10829b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 10839b94acceSBarry Smith @*/ 10848ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 10859b94acceSBarry Smith { 1086272ac6f2SLois Curfman McInnes int ierr, flg; 108777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 108877c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 10898ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 10909b94acceSBarry Smith 1091c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1092c1f60f51SBarry Smith /* 1093c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1094dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1095c1f60f51SBarry Smith */ 1096c1f60f51SBarry Smith if (flg) { 1097c1f60f51SBarry Smith Mat J; 1098c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1099c1f60f51SBarry Smith PLogObjectParent(snes,J); 1100c1f60f51SBarry Smith snes->mfshell = J; 1101c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1102c1f60f51SBarry Smith snes->jacobian = J; 1103c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1104c1f60f51SBarry Smith } 1105c1f60f51SBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1106c1f60f51SBarry Smith snes->jacobian = J; 1107c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1108e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free operator option"); 1109c1f60f51SBarry Smith } 1110272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1111c1f60f51SBarry Smith /* 1112dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1113c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1114c1f60f51SBarry Smith */ 111531615d04SLois Curfman McInnes if (flg) { 1116272ac6f2SLois Curfman McInnes Mat J; 1117272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1118272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1119272ac6f2SLois Curfman McInnes snes->mfshell = J; 112083e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 112183e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 112294a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 112383e56afdSLois Curfman McInnes } 112483e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 112583e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 112694a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1127e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free option"); 1128272ac6f2SLois Curfman McInnes } 11299b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1130e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1131e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1132e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first"); 1133e3372554SBarry Smith if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector"); 1134a703fe33SLois Curfman McInnes 1135a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 113640191667SLois Curfman McInnes if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1137a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1138a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1139a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1140a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1141a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1142a703fe33SLois Curfman McInnes } 11439b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1144e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first"); 1145e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first"); 114637fcc0dbSBarry Smith if (!snes->computeumfunction) 1147e3372554SBarry Smith SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first"); 1148e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first"); 1149e3372554SBarry Smith } else SETERRQ(1,0,"Unknown method class"); 1150a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1151a703fe33SLois Curfman McInnes snes->setup_called = 1; 1152a703fe33SLois Curfman McInnes return 0; 11539b94acceSBarry Smith } 11549b94acceSBarry Smith 11555615d1e5SSatish Balay #undef __FUNC__ 11565615d1e5SSatish Balay #define __FUNC__ "SNESDestroy" 11579b94acceSBarry Smith /*@C 11589b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 11599b94acceSBarry Smith with SNESCreate(). 11609b94acceSBarry Smith 11619b94acceSBarry Smith Input Parameter: 11629b94acceSBarry Smith . snes - the SNES context 11639b94acceSBarry Smith 11649b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 11659b94acceSBarry Smith 116663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 11679b94acceSBarry Smith @*/ 11689b94acceSBarry Smith int SNESDestroy(SNES snes) 11699b94acceSBarry Smith { 11709b94acceSBarry Smith int ierr; 117177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 11729b94acceSBarry Smith ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr); 11730452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 11749b94acceSBarry Smith if (snes->mfshell) MatDestroy(snes->mfshell); 11759b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 11763e2e452bSBarry Smith if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 11776f24a144SLois Curfman McInnes if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 11789b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 11790452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 11809b94acceSBarry Smith return 0; 11819b94acceSBarry Smith } 11829b94acceSBarry Smith 11839b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 11849b94acceSBarry Smith 11855615d1e5SSatish Balay #undef __FUNC__ 11865615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 11879b94acceSBarry Smith /*@ 1188d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 11899b94acceSBarry Smith 11909b94acceSBarry Smith Input Parameters: 11919b94acceSBarry Smith . snes - the SNES context 119233174efeSLois Curfman McInnes . atol - absolute convergence tolerance 119333174efeSLois Curfman McInnes . rtol - relative convergence tolerance 119433174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 119533174efeSLois Curfman McInnes of the change in the solution between steps 119633174efeSLois Curfman McInnes . maxit - maximum number of iterations 119733174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 11989b94acceSBarry Smith 119933174efeSLois Curfman McInnes Options Database Keys: 120033174efeSLois Curfman McInnes $ -snes_atol <atol> 120133174efeSLois Curfman McInnes $ -snes_rtol <rtol> 120233174efeSLois Curfman McInnes $ -snes_stol <stol> 120333174efeSLois Curfman McInnes $ -snes_max_it <maxit> 120433174efeSLois Curfman McInnes $ -snes_max_funcs <maxf> 12059b94acceSBarry Smith 1206d7a720efSLois Curfman McInnes Notes: 12079b94acceSBarry Smith The default maximum number of iterations is 50. 12089b94acceSBarry Smith The default maximum number of function evaluations is 1000. 12099b94acceSBarry Smith 121033174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 12119b94acceSBarry Smith 1212d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 12139b94acceSBarry Smith @*/ 1214d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 12159b94acceSBarry Smith { 121677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1217d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1218d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1219d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1220d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1221d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 12229b94acceSBarry Smith return 0; 12239b94acceSBarry Smith } 12249b94acceSBarry Smith 12255615d1e5SSatish Balay #undef __FUNC__ 12265615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 12279b94acceSBarry Smith /*@ 122833174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 122933174efeSLois Curfman McInnes 123033174efeSLois Curfman McInnes Input Parameters: 123133174efeSLois Curfman McInnes . snes - the SNES context 123233174efeSLois Curfman McInnes . atol - absolute convergence tolerance 123333174efeSLois Curfman McInnes . rtol - relative convergence tolerance 123433174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 123533174efeSLois Curfman McInnes of the change in the solution between steps 123633174efeSLois Curfman McInnes . maxit - maximum number of iterations 123733174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 123833174efeSLois Curfman McInnes 123933174efeSLois Curfman McInnes Notes: 124033174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 124133174efeSLois Curfman McInnes 124233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 124333174efeSLois Curfman McInnes 124433174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 124533174efeSLois Curfman McInnes @*/ 124633174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 124733174efeSLois Curfman McInnes { 124833174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 124933174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 125033174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 125133174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 125233174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 125333174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 125433174efeSLois Curfman McInnes return 0; 125533174efeSLois Curfman McInnes } 125633174efeSLois Curfman McInnes 12575615d1e5SSatish Balay #undef __FUNC__ 12585615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 125933174efeSLois Curfman McInnes /*@ 12609b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 12619b94acceSBarry Smith 12629b94acceSBarry Smith Input Parameters: 12639b94acceSBarry Smith . snes - the SNES context 12649b94acceSBarry Smith . tol - tolerance 12659b94acceSBarry Smith 12669b94acceSBarry Smith Options Database Key: 1267d7a720efSLois Curfman McInnes $ -snes_trtol <tol> 12689b94acceSBarry Smith 12699b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 12709b94acceSBarry Smith 1271d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 12729b94acceSBarry Smith @*/ 12739b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 12749b94acceSBarry Smith { 127577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12769b94acceSBarry Smith snes->deltatol = tol; 12779b94acceSBarry Smith return 0; 12789b94acceSBarry Smith } 12799b94acceSBarry Smith 12805615d1e5SSatish Balay #undef __FUNC__ 12815615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 12829b94acceSBarry Smith /*@ 128377c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 12849b94acceSBarry Smith for unconstrained minimization solvers. 12859b94acceSBarry Smith 12869b94acceSBarry Smith Input Parameters: 12879b94acceSBarry Smith . snes - the SNES context 12889b94acceSBarry Smith . ftol - minimum function tolerance 12899b94acceSBarry Smith 12909b94acceSBarry Smith Options Database Key: 1291d7a720efSLois Curfman McInnes $ -snes_fmin <ftol> 12929b94acceSBarry Smith 12939b94acceSBarry Smith Note: 129477c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 12959b94acceSBarry Smith methods only. 12969b94acceSBarry Smith 12979b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 12989b94acceSBarry Smith 1299d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 13009b94acceSBarry Smith @*/ 130177c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 13029b94acceSBarry Smith { 130377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13049b94acceSBarry Smith snes->fmin = ftol; 13059b94acceSBarry Smith return 0; 13069b94acceSBarry Smith } 13079b94acceSBarry Smith 13089b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 13099b94acceSBarry Smith 13105615d1e5SSatish Balay #undef __FUNC__ 13115615d1e5SSatish Balay #define __FUNC__ "SNESSetMonitor" 13129b94acceSBarry Smith /*@C 13139b94acceSBarry Smith SNESSetMonitor - Sets the function that is to be used at every 13149b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 13159b94acceSBarry Smith progress. 13169b94acceSBarry Smith 13179b94acceSBarry Smith Input Parameters: 13189b94acceSBarry Smith . snes - the SNES context 13199b94acceSBarry Smith . func - monitoring routine 1320044dda88SLois Curfman McInnes . mctx - [optional] user-defined context for private data for the 1321044dda88SLois Curfman McInnes monitor routine (may be PETSC_NULL) 13229b94acceSBarry Smith 13239b94acceSBarry Smith Calling sequence of func: 1324682d7d0cSBarry Smith int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 13259b94acceSBarry Smith 13269b94acceSBarry Smith $ snes - the SNES context 13279b94acceSBarry Smith $ its - iteration number 1328044dda88SLois Curfman McInnes $ mctx - [optional] monitoring context 13299b94acceSBarry Smith $ 13309b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13319b94acceSBarry Smith $ norm - 2-norm function value (may be estimated) 13329b94acceSBarry Smith $ 13339b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13349b94acceSBarry Smith $ norm - 2-norm gradient value (may be estimated) 13359b94acceSBarry Smith 13369665c990SLois Curfman McInnes Options Database Keys: 13379665c990SLois Curfman McInnes $ -snes_monitor : sets SNESDefaultMonitor() 13389665c990SLois Curfman McInnes $ -snes_xmonitor : sets line graph monitor, 13399665c990SLois Curfman McInnes $ uses SNESLGMonitorCreate() 13409665c990SLois Curfman McInnes $ -snes_cancelmonitors : cancels all monitors that have 13419665c990SLois Curfman McInnes $ been hardwired into a code by 13429665c990SLois Curfman McInnes $ calls to SNESSetMonitor(), but 13439665c990SLois Curfman McInnes $ does not cancel those set via 13449665c990SLois Curfman McInnes $ the options database. 13459665c990SLois Curfman McInnes 13469665c990SLois Curfman McInnes 1347639f9d9dSBarry Smith Notes: 13486bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 13496bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 13506bc08f3fSLois Curfman McInnes order in which they were set. 1351639f9d9dSBarry Smith 13529b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 13539b94acceSBarry Smith 13549b94acceSBarry Smith .seealso: SNESDefaultMonitor() 13559b94acceSBarry Smith @*/ 135674679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 13579b94acceSBarry Smith { 1358639f9d9dSBarry Smith if (!func) { 1359639f9d9dSBarry Smith snes->numbermonitors = 0; 1360639f9d9dSBarry Smith return 0; 1361639f9d9dSBarry Smith } 1362639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1363e3372554SBarry Smith SETERRQ(1,0,"Too many monitors set"); 1364639f9d9dSBarry Smith } 1365639f9d9dSBarry Smith 1366639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1367639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 13689b94acceSBarry Smith return 0; 13699b94acceSBarry Smith } 13709b94acceSBarry Smith 13715615d1e5SSatish Balay #undef __FUNC__ 13725615d1e5SSatish Balay #define __FUNC__ "SNESSetConvergenceTest" 13739b94acceSBarry Smith /*@C 13749b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13759b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13769b94acceSBarry Smith 13779b94acceSBarry Smith Input Parameters: 13789b94acceSBarry Smith . snes - the SNES context 13799b94acceSBarry Smith . func - routine to test for convergence 1380044dda88SLois Curfman McInnes . cctx - [optional] context for private data for the convergence routine 1381044dda88SLois Curfman McInnes (may be PETSC_NULL) 13829b94acceSBarry Smith 13839b94acceSBarry Smith Calling sequence of func: 13849b94acceSBarry Smith int func (SNES snes,double xnorm,double gnorm, 13859b94acceSBarry Smith double f,void *cctx) 13869b94acceSBarry Smith 13879b94acceSBarry Smith $ snes - the SNES context 1388044dda88SLois Curfman McInnes $ cctx - [optional] convergence context 13899b94acceSBarry Smith $ xnorm - 2-norm of current iterate 13909b94acceSBarry Smith $ 13919b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13929b94acceSBarry Smith $ gnorm - 2-norm of current step 13939b94acceSBarry Smith $ f - 2-norm of function 13949b94acceSBarry Smith $ 13959b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13969b94acceSBarry Smith $ gnorm - 2-norm of current gradient 13979b94acceSBarry Smith $ f - function value 13989b94acceSBarry Smith 13999b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 14009b94acceSBarry Smith 140140191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 140240191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 14039b94acceSBarry Smith @*/ 140474679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 14059b94acceSBarry Smith { 14069b94acceSBarry Smith (snes)->converged = func; 14079b94acceSBarry Smith (snes)->cnvP = cctx; 14089b94acceSBarry Smith return 0; 14099b94acceSBarry Smith } 14109b94acceSBarry Smith 14115615d1e5SSatish Balay #undef __FUNC__ 14125615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 14139b94acceSBarry Smith /* 14149b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 14159b94acceSBarry Smith positive parameter delta. 14169b94acceSBarry Smith 14179b94acceSBarry Smith Input Parameters: 14189b94acceSBarry Smith . snes - the SNES context 14199b94acceSBarry Smith . y - approximate solution of linear system 14209b94acceSBarry Smith . fnorm - 2-norm of current function 14219b94acceSBarry Smith . delta - trust region size 14229b94acceSBarry Smith 14239b94acceSBarry Smith Output Parameters: 14249b94acceSBarry Smith . gpnorm - predicted function norm at the new point, assuming local 14259b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 14269b94acceSBarry Smith region, and exceeds zero otherwise. 14279b94acceSBarry Smith . ynorm - 2-norm of the step 14289b94acceSBarry Smith 14299b94acceSBarry Smith Note: 143040191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 14319b94acceSBarry Smith is set to be the maximum allowable step size. 14329b94acceSBarry Smith 14339b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 14349b94acceSBarry Smith */ 14359b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 14369b94acceSBarry Smith double *gpnorm,double *ynorm) 14379b94acceSBarry Smith { 14389b94acceSBarry Smith double norm; 14399b94acceSBarry Smith Scalar cnorm; 1440cddf8d76SBarry Smith VecNorm(y,NORM_2, &norm ); 14419b94acceSBarry Smith if (norm > *delta) { 14429b94acceSBarry Smith norm = *delta/norm; 14439b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 14449b94acceSBarry Smith cnorm = norm; 14459b94acceSBarry Smith VecScale( &cnorm, y ); 14469b94acceSBarry Smith *ynorm = *delta; 14479b94acceSBarry Smith } else { 14489b94acceSBarry Smith *gpnorm = 0.0; 14499b94acceSBarry Smith *ynorm = norm; 14509b94acceSBarry Smith } 14519b94acceSBarry Smith return 0; 14529b94acceSBarry Smith } 14539b94acceSBarry Smith 14545615d1e5SSatish Balay #undef __FUNC__ 14555615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 14569b94acceSBarry Smith /*@ 14579b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 145863a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 14599b94acceSBarry Smith 14609b94acceSBarry Smith Input Parameter: 14619b94acceSBarry Smith . snes - the SNES context 14628ddd3da0SLois Curfman McInnes . x - the solution vector 14639b94acceSBarry Smith 14649b94acceSBarry Smith Output Parameter: 14659b94acceSBarry Smith its - number of iterations until termination 14669b94acceSBarry Smith 14678ddd3da0SLois Curfman McInnes Note: 14688ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 14698ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 14708ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 14718ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 14728ddd3da0SLois Curfman McInnes 14739b94acceSBarry Smith .keywords: SNES, nonlinear, solve 14749b94acceSBarry Smith 147563a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 14769b94acceSBarry Smith @*/ 14778ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 14789b94acceSBarry Smith { 14793c7409f5SSatish Balay int ierr, flg; 1480052efed2SBarry Smith 148177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 148274679c65SBarry Smith PetscValidIntPointer(its); 1483c4fc05e7SBarry Smith if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1484c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 14859b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 14869b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 14879b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 14883c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 14896d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 14909b94acceSBarry Smith return 0; 14919b94acceSBarry Smith } 14929b94acceSBarry Smith 14939b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 149437fcc0dbSBarry Smith static NRList *__SNESList = 0; 14959b94acceSBarry Smith 14965615d1e5SSatish Balay #undef __FUNC__ 14975615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 14989b94acceSBarry Smith /*@ 14994b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 15009b94acceSBarry Smith 15019b94acceSBarry Smith Input Parameters: 15029b94acceSBarry Smith . snes - the SNES context 15039b94acceSBarry Smith . method - a known method 15049b94acceSBarry Smith 1505ae12b187SLois Curfman McInnes Options Database Command: 1506ae12b187SLois Curfman McInnes $ -snes_type <method> 1507ae12b187SLois Curfman McInnes $ Use -help for a list of available methods 1508ae12b187SLois Curfman McInnes $ (for instance, ls or tr) 1509ae12b187SLois Curfman McInnes 15109b94acceSBarry Smith Notes: 15119b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 15129b94acceSBarry Smith $ Systems of nonlinear equations: 151340191667SLois Curfman McInnes $ SNES_EQ_LS - Newton's method with line search 151440191667SLois Curfman McInnes $ SNES_EQ_TR - Newton's method with trust region 15159b94acceSBarry Smith $ Unconstrained minimization: 151640191667SLois Curfman McInnes $ SNES_UM_TR - Newton's method with trust region 151740191667SLois Curfman McInnes $ SNES_UM_LS - Newton's method with line search 15189b94acceSBarry Smith 1519ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1520ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1521ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1522ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1523ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1524ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1525ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1526ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1527ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1528ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1529a703fe33SLois Curfman McInnes 1530f536c699SSatish Balay .keywords: SNES, set, method 15319b94acceSBarry Smith @*/ 15324b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 15339b94acceSBarry Smith { 15349b94acceSBarry Smith int (*r)(SNES); 153577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 15367d1a2b2bSBarry Smith if (snes->type == (int) method) return 0; 15377d1a2b2bSBarry Smith 15389b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 153937fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 1540e3372554SBarry Smith if (!__SNESList) {SETERRQ(1,0,"Could not get methods");} 154137fcc0dbSBarry Smith r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1542e3372554SBarry Smith if (!r) {SETERRQ(1,0,"Unknown method");} 15430452661fSBarry Smith if (snes->data) PetscFree(snes->data); 15449b94acceSBarry Smith snes->set_method_called = 1; 15459b94acceSBarry Smith return (*r)(snes); 15469b94acceSBarry Smith } 15479b94acceSBarry Smith 15489b94acceSBarry Smith /* --------------------------------------------------------------------- */ 15495615d1e5SSatish Balay #undef __FUNC__ 15505615d1e5SSatish Balay #define __FUNC__ "SNESRegister" 15519b94acceSBarry Smith /*@C 15529b94acceSBarry Smith SNESRegister - Adds the method to the nonlinear solver package, given 15534b0e389bSBarry Smith a function pointer and a nonlinear solver name of the type SNESType. 15549b94acceSBarry Smith 15559b94acceSBarry Smith Input Parameters: 155640191667SLois Curfman McInnes . name - for instance SNES_EQ_LS, SNES_EQ_TR, ... 155740191667SLois Curfman McInnes . sname - corresponding string for name 15589b94acceSBarry Smith . create - routine to create method context 15599b94acceSBarry Smith 15609b94acceSBarry Smith .keywords: SNES, nonlinear, register 15619b94acceSBarry Smith 15629b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy() 15639b94acceSBarry Smith @*/ 15649b94acceSBarry Smith int SNESRegister(int name, char *sname, int (*create)(SNES)) 15659b94acceSBarry Smith { 15669b94acceSBarry Smith int ierr; 156737fcc0dbSBarry Smith if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 156837fcc0dbSBarry Smith NRRegister( __SNESList, name, sname, (int (*)(void*))create ); 15699b94acceSBarry Smith return 0; 15709b94acceSBarry Smith } 1571a847f771SSatish Balay 15729b94acceSBarry Smith /* --------------------------------------------------------------------- */ 15735615d1e5SSatish Balay #undef __FUNC__ 15745615d1e5SSatish Balay #define __FUNC__ "SNESRegisterDestroy" 15759b94acceSBarry Smith /*@C 15769b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 15779b94acceSBarry Smith registered by SNESRegister(). 15789b94acceSBarry Smith 15799b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 15809b94acceSBarry Smith 15819b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 15829b94acceSBarry Smith @*/ 15839b94acceSBarry Smith int SNESRegisterDestroy() 15849b94acceSBarry Smith { 158537fcc0dbSBarry Smith if (__SNESList) { 158637fcc0dbSBarry Smith NRDestroy( __SNESList ); 158737fcc0dbSBarry Smith __SNESList = 0; 15889b94acceSBarry Smith } 15899b94acceSBarry Smith return 0; 15909b94acceSBarry Smith } 15919b94acceSBarry Smith 15925615d1e5SSatish Balay #undef __FUNC__ 15935615d1e5SSatish Balay #define __FUNC__ "SNESGetTypeFromOptions_Private" 15949b94acceSBarry Smith /* 15954b0e389bSBarry Smith SNESGetTypeFromOptions_Private - Sets the selected method from the 15969b94acceSBarry Smith options database. 15979b94acceSBarry Smith 15989b94acceSBarry Smith Input Parameter: 15999b94acceSBarry Smith . ctx - the SNES context 16009b94acceSBarry Smith 16019b94acceSBarry Smith Output Parameter: 16029b94acceSBarry Smith . method - solver method 16039b94acceSBarry Smith 16049b94acceSBarry Smith Returns: 16059b94acceSBarry Smith Returns 1 if the method is found; 0 otherwise. 16069b94acceSBarry Smith 16079b94acceSBarry Smith Options Database Key: 16084b0e389bSBarry Smith $ -snes_type method 16099b94acceSBarry Smith */ 1610052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 16119b94acceSBarry Smith { 1612052efed2SBarry Smith int ierr; 16139b94acceSBarry Smith char sbuf[50]; 16145baf8537SBarry Smith 1615052efed2SBarry Smith ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1616052efed2SBarry Smith if (*flg) { 1617052efed2SBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 161837fcc0dbSBarry Smith *method = (SNESType)NRFindID( __SNESList, sbuf ); 16199b94acceSBarry Smith } 16209b94acceSBarry Smith return 0; 16219b94acceSBarry Smith } 16229b94acceSBarry Smith 16235615d1e5SSatish Balay #undef __FUNC__ 16245615d1e5SSatish Balay #define __FUNC__ "SNESGetType" 16259b94acceSBarry Smith /*@C 16269a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 16279b94acceSBarry Smith 16289b94acceSBarry Smith Input Parameter: 16294b0e389bSBarry Smith . snes - nonlinear solver context 16309b94acceSBarry Smith 16319b94acceSBarry Smith Output Parameter: 16329a28b0a6SLois Curfman McInnes . method - SNES method (or use PETSC_NULL) 16339a28b0a6SLois Curfman McInnes . name - name of SNES method (or use PETSC_NULL) 16349b94acceSBarry Smith 16359b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 16369b94acceSBarry Smith @*/ 16374b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name) 16389b94acceSBarry Smith { 163906a9b0adSLois Curfman McInnes int ierr; 164037fcc0dbSBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 16414b0e389bSBarry Smith if (method) *method = (SNESType) snes->type; 164237fcc0dbSBarry Smith if (name) *name = NRFindName( __SNESList, (int) snes->type ); 16439b94acceSBarry Smith return 0; 16449b94acceSBarry Smith } 16459b94acceSBarry Smith 16469b94acceSBarry Smith #include <stdio.h> 16475615d1e5SSatish Balay #undef __FUNC__ 16485615d1e5SSatish Balay #define __FUNC__ "SNESPrintTypes_Private" 16499b94acceSBarry Smith /* 16504b0e389bSBarry Smith SNESPrintTypes_Private - Prints the SNES methods available from the 16519b94acceSBarry Smith options database. 16529b94acceSBarry Smith 16539b94acceSBarry Smith Input Parameters: 165433455573SSatish Balay . comm - communicator (usually MPI_COMM_WORLD) 16559b94acceSBarry Smith . prefix - prefix (usually "-") 16564b0e389bSBarry Smith . name - the options database name (by default "snes_type") 16579b94acceSBarry Smith */ 165833455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 16599b94acceSBarry Smith { 16609b94acceSBarry Smith FuncList *entry; 166137fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 166237fcc0dbSBarry Smith entry = __SNESList->head; 166377c4ece6SBarry Smith PetscPrintf(comm," %s%s (one of)",prefix,name); 16649b94acceSBarry Smith while (entry) { 166577c4ece6SBarry Smith PetscPrintf(comm," %s",entry->name); 16669b94acceSBarry Smith entry = entry->next; 16679b94acceSBarry Smith } 166877c4ece6SBarry Smith PetscPrintf(comm,"\n"); 16699b94acceSBarry Smith return 0; 16709b94acceSBarry Smith } 16719b94acceSBarry Smith 16725615d1e5SSatish Balay #undef __FUNC__ 16735615d1e5SSatish Balay #define __FUNC__ "SNESGetSolution" 16749b94acceSBarry Smith /*@C 16759b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 16769b94acceSBarry Smith stored. 16779b94acceSBarry Smith 16789b94acceSBarry Smith Input Parameter: 16799b94acceSBarry Smith . snes - the SNES context 16809b94acceSBarry Smith 16819b94acceSBarry Smith Output Parameter: 16829b94acceSBarry Smith . x - the solution 16839b94acceSBarry Smith 16849b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 16859b94acceSBarry Smith 1686a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 16879b94acceSBarry Smith @*/ 16889b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 16899b94acceSBarry Smith { 169077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16919b94acceSBarry Smith *x = snes->vec_sol_always; 16929b94acceSBarry Smith return 0; 16939b94acceSBarry Smith } 16949b94acceSBarry Smith 16955615d1e5SSatish Balay #undef __FUNC__ 16965615d1e5SSatish Balay #define __FUNC__ "SNESGetSolutionUpdate" 16979b94acceSBarry Smith /*@C 16989b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 16999b94acceSBarry Smith stored. 17009b94acceSBarry Smith 17019b94acceSBarry Smith Input Parameter: 17029b94acceSBarry Smith . snes - the SNES context 17039b94acceSBarry Smith 17049b94acceSBarry Smith Output Parameter: 17059b94acceSBarry Smith . x - the solution update 17069b94acceSBarry Smith 17079b94acceSBarry Smith Notes: 17089b94acceSBarry Smith This vector is implementation dependent. 17099b94acceSBarry Smith 17109b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 17119b94acceSBarry Smith 17129b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 17139b94acceSBarry Smith @*/ 17149b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 17159b94acceSBarry Smith { 171677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17179b94acceSBarry Smith *x = snes->vec_sol_update_always; 17189b94acceSBarry Smith return 0; 17199b94acceSBarry Smith } 17209b94acceSBarry Smith 17215615d1e5SSatish Balay #undef __FUNC__ 17225615d1e5SSatish Balay #define __FUNC__ "SNESGetFunction" 17239b94acceSBarry Smith /*@C 17243638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 17259b94acceSBarry Smith 17269b94acceSBarry Smith Input Parameter: 17279b94acceSBarry Smith . snes - the SNES context 17289b94acceSBarry Smith 17299b94acceSBarry Smith Output Parameter: 17303638b69dSLois Curfman McInnes . r - the function 17319b94acceSBarry Smith 17329b94acceSBarry Smith Notes: 17339b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 17349b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 17359b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 17369b94acceSBarry Smith 1737a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 17389b94acceSBarry Smith 173931615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 174031615d04SLois Curfman McInnes SNESGetGradient() 17419b94acceSBarry Smith @*/ 17429b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 17439b94acceSBarry Smith { 174477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1745e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0, 17467c792091SSatish Balay "For SNES_NONLINEAR_EQUATIONS only"); 17479b94acceSBarry Smith *r = snes->vec_func_always; 17489b94acceSBarry Smith return 0; 17499b94acceSBarry Smith } 17509b94acceSBarry Smith 17515615d1e5SSatish Balay #undef __FUNC__ 17525615d1e5SSatish Balay #define __FUNC__ "SNESGetGradient" 17539b94acceSBarry Smith /*@C 17543638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 17559b94acceSBarry Smith 17569b94acceSBarry Smith Input Parameter: 17579b94acceSBarry Smith . snes - the SNES context 17589b94acceSBarry Smith 17599b94acceSBarry Smith Output Parameter: 17609b94acceSBarry Smith . r - the gradient 17619b94acceSBarry Smith 17629b94acceSBarry Smith Notes: 17639b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 17649b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 17659b94acceSBarry Smith SNESGetFunction(). 17669b94acceSBarry Smith 17679b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 17689b94acceSBarry Smith 176931615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 17709b94acceSBarry Smith @*/ 17719b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 17729b94acceSBarry Smith { 177377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1774e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 177563c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 17769b94acceSBarry Smith *r = snes->vec_func_always; 17779b94acceSBarry Smith return 0; 17789b94acceSBarry Smith } 17799b94acceSBarry Smith 17805615d1e5SSatish Balay #undef __FUNC__ 17815615d1e5SSatish Balay #define __FUNC__ "SNESGetMinimizationFunction" 17829b94acceSBarry Smith /*@ 17839b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 17849b94acceSBarry Smith unconstrained minimization problems. 17859b94acceSBarry Smith 17869b94acceSBarry Smith Input Parameter: 17879b94acceSBarry Smith . snes - the SNES context 17889b94acceSBarry Smith 17899b94acceSBarry Smith Output Parameter: 17909b94acceSBarry Smith . r - the function 17919b94acceSBarry Smith 17929b94acceSBarry Smith Notes: 17939b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 17949b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 17959b94acceSBarry Smith SNESGetFunction(). 17969b94acceSBarry Smith 17979b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 17989b94acceSBarry Smith 179931615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 18009b94acceSBarry Smith @*/ 18019b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 18029b94acceSBarry Smith { 180377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 180474679c65SBarry Smith PetscValidScalarPointer(r); 1805e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 180663c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 18079b94acceSBarry Smith *r = snes->fc; 18089b94acceSBarry Smith return 0; 18099b94acceSBarry Smith } 18109b94acceSBarry Smith 18115615d1e5SSatish Balay #undef __FUNC__ 18125615d1e5SSatish Balay #define __FUNC__ "SNESSetOptionsPrefix" 18133c7409f5SSatish Balay /*@C 18143c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1815639f9d9dSBarry Smith SNES options in the database. You must NOT include the - at the beginning of 1816639f9d9dSBarry Smith the prefix name. 18173c7409f5SSatish Balay 18183c7409f5SSatish Balay Input Parameter: 18193c7409f5SSatish Balay . snes - the SNES context 18203c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 18213c7409f5SSatish Balay 18223c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1823a86d99e1SLois Curfman McInnes 1824a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 18253c7409f5SSatish Balay @*/ 18263c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 18273c7409f5SSatish Balay { 18283c7409f5SSatish Balay int ierr; 18293c7409f5SSatish Balay 183077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1831639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 18323c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 18333c7409f5SSatish Balay return 0; 18343c7409f5SSatish Balay } 18353c7409f5SSatish Balay 18365615d1e5SSatish Balay #undef __FUNC__ 18375615d1e5SSatish Balay #define __FUNC__ "SNESAppendOptionsPrefix" 18383c7409f5SSatish Balay /*@C 1839f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1840639f9d9dSBarry Smith SNES options in the database. You must NOT include the - at the beginning of 1841639f9d9dSBarry Smith the prefix name. 18423c7409f5SSatish Balay 18433c7409f5SSatish Balay Input Parameter: 18443c7409f5SSatish Balay . snes - the SNES context 18453c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 18463c7409f5SSatish Balay 18473c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1848a86d99e1SLois Curfman McInnes 1849a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 18503c7409f5SSatish Balay @*/ 18513c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 18523c7409f5SSatish Balay { 18533c7409f5SSatish Balay int ierr; 18543c7409f5SSatish Balay 185577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1856639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 18573c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 18583c7409f5SSatish Balay return 0; 18593c7409f5SSatish Balay } 18603c7409f5SSatish Balay 18615615d1e5SSatish Balay #undef __FUNC__ 18625615d1e5SSatish Balay #define __FUNC__ "SNESGetOptionsPrefix" 1863c83590e2SSatish Balay /*@ 18643c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 18653c7409f5SSatish Balay SNES options in the database. 18663c7409f5SSatish Balay 18673c7409f5SSatish Balay Input Parameter: 18683c7409f5SSatish Balay . snes - the SNES context 18693c7409f5SSatish Balay 18703c7409f5SSatish Balay Output Parameter: 18713c7409f5SSatish Balay . prefix - pointer to the prefix string used 18723c7409f5SSatish Balay 18733c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 1874a86d99e1SLois Curfman McInnes 1875a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 18763c7409f5SSatish Balay @*/ 18773c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 18783c7409f5SSatish Balay { 18793c7409f5SSatish Balay int ierr; 18803c7409f5SSatish Balay 188177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1882639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 18833c7409f5SSatish Balay return 0; 18843c7409f5SSatish Balay } 18853c7409f5SSatish Balay 18863c7409f5SSatish Balay 18879b94acceSBarry Smith 18889b94acceSBarry Smith 18899b94acceSBarry Smith 1890