19b94acceSBarry Smith #ifndef lint 2*81f57ec7SBarry Smith static char vcid[] = "$Id: snes.c,v 1.99 1996/11/13 16:25:34 curfman Exp bsmith $"; 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 149b94acceSBarry Smith /*@ 159b94acceSBarry Smith SNESView - Prints the SNES data structure. 169b94acceSBarry Smith 179b94acceSBarry Smith Input Parameters: 189b94acceSBarry Smith . SNES - the SNES context 199b94acceSBarry Smith . viewer - visualization context 209b94acceSBarry Smith 219b94acceSBarry Smith Options Database Key: 229b94acceSBarry Smith $ -snes_view : calls SNESView() at end of SNESSolve() 239b94acceSBarry Smith 249b94acceSBarry Smith Notes: 259b94acceSBarry Smith The available visualization contexts include 266d4a8577SBarry Smith $ VIEWER_STDOUT_SELF - standard output (default) 276d4a8577SBarry Smith $ VIEWER_STDOUT_WORLD - synchronized standard 289b94acceSBarry Smith $ output where only the first processor opens 299b94acceSBarry Smith $ the file. All other processors send their 309b94acceSBarry Smith $ data to the first processor to print. 319b94acceSBarry Smith 329b94acceSBarry Smith The user can open alternative vistualization contexts with 33dbb450caSBarry Smith $ ViewerFileOpenASCII() - output to a specified file 349b94acceSBarry Smith 359b94acceSBarry Smith .keywords: SNES, view 369b94acceSBarry Smith 37dbb450caSBarry Smith .seealso: ViewerFileOpenASCII() 389b94acceSBarry Smith @*/ 399b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 409b94acceSBarry Smith { 419b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 429b94acceSBarry Smith FILE *fd; 439b94acceSBarry Smith int ierr; 449b94acceSBarry Smith SLES sles; 459b94acceSBarry Smith char *method; 4619bcc07fSBarry Smith ViewerType vtype; 479b94acceSBarry Smith 4874679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 4974679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5074679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5174679c65SBarry Smith 5219bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 5319bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 5490ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 5577c4ece6SBarry Smith PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 564b0e389bSBarry Smith SNESGetType(snes,PETSC_NULL,&method); 5777c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," method: %s\n",method); 589b94acceSBarry Smith if (snes->view) (*snes->view)((PetscObject)snes,viewer); 5977c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 609b94acceSBarry Smith " maximum iterations=%d, maximum function evaluations=%d\n", 619b94acceSBarry Smith snes->max_its,snes->max_funcs); 6277c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 639b94acceSBarry Smith " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 649b94acceSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol); 659b94acceSBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 6677c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 679b94acceSBarry Smith if (snes->ksp_ewconv) { 689b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 699b94acceSBarry Smith if (kctx) { 7077c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 719b94acceSBarry Smith " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 729b94acceSBarry Smith kctx->version); 7377c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 749b94acceSBarry Smith " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 759b94acceSBarry Smith kctx->rtol_max,kctx->threshold); 7677c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 779b94acceSBarry Smith kctx->gamma,kctx->alpha,kctx->alpha2); 789b94acceSBarry Smith } 799b94acceSBarry Smith } 80c30f285eSLois Curfman McInnes } else if (vtype == STRING_VIEWER) { 81c30f285eSLois Curfman McInnes SNESGetType(snes,PETSC_NULL,&method); 82c30f285eSLois Curfman McInnes ViewerStringSPrintf(viewer," %-3.3s",method); 8319bcc07fSBarry Smith } 849b94acceSBarry Smith SNESGetSLES(snes,&sles); 859b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 869b94acceSBarry Smith return 0; 879b94acceSBarry Smith } 889b94acceSBarry Smith 89639f9d9dSBarry Smith /* 90639f9d9dSBarry Smith We retain a list of functions that also take SNES command 91639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 92639f9d9dSBarry Smith */ 93639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 94639f9d9dSBarry Smith static int numberofsetfromoptions; 95639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 96639f9d9dSBarry Smith 97639f9d9dSBarry Smith /*@ 98639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 99639f9d9dSBarry Smith 100639f9d9dSBarry Smith Input Parameter: 101639f9d9dSBarry Smith . snescheck - function that checks for options 102639f9d9dSBarry Smith 103639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 104639f9d9dSBarry Smith @*/ 105639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 106639f9d9dSBarry Smith { 107639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 108639f9d9dSBarry Smith SETERRQ(1,"SNESAddOptionsChecker:Too many options checkers, only 5 allowed"); 109639f9d9dSBarry Smith } 110639f9d9dSBarry Smith 111639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 112639f9d9dSBarry Smith return 0; 113639f9d9dSBarry Smith } 114639f9d9dSBarry Smith 1159b94acceSBarry Smith /*@ 116682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1179b94acceSBarry Smith 1189b94acceSBarry Smith Input Parameter: 1199b94acceSBarry Smith . snes - the SNES context 1209b94acceSBarry Smith 1219b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1229b94acceSBarry Smith 123a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1249b94acceSBarry Smith @*/ 1259b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1269b94acceSBarry Smith { 1274b0e389bSBarry Smith SNESType method; 1289b94acceSBarry Smith double tmp; 1299b94acceSBarry Smith SLES sles; 130*81f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1313c7409f5SSatish Balay int version = PETSC_DEFAULT; 1329b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1339b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1349b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1359b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1369b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1379b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1389b94acceSBarry Smith 139*81f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 140*81f57ec7SBarry Smith 141*81f57ec7SBarry Smith 14277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 14374679c65SBarry Smith if (snes->setup_called) SETERRQ(1,"SNESSetFromOptions:Must call prior to SNESSetUp"); 144052efed2SBarry Smith ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 145052efed2SBarry Smith if (flg) { 146052efed2SBarry Smith ierr = SNESSetType(snes,method); CHKERRQ(ierr); 1479b94acceSBarry Smith } 1485334005bSBarry Smith else if (!snes->set_method_called) { 1495334005bSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 15040191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 1515334005bSBarry Smith } 1525334005bSBarry Smith else { 15340191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 1545334005bSBarry Smith } 1555334005bSBarry Smith } 1565334005bSBarry Smith 1573c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1583c7409f5SSatish Balay if (flg) { SNESPrintHelp(snes); } 1593c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 160d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);} 1613c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 162d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);} 1633c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 164d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); } 1653c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 1663c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 167d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 168d7a720efSLois Curfman McInnes if (flg) { SNESSetTrustRegionTolerance(snes,tmp); } 169d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 170d7a720efSLois Curfman McInnes if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); } 1713c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 1723c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 173b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 174b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 175b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 176b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 177b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 178b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 179b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 1809b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 1819b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 1825bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 1835bbad29bSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);} 1843c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 185639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 1863c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 187d31a3109SLois Curfman McInnes if (flg) { 188639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 189d31a3109SLois Curfman McInnes } 190*81f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 1913c7409f5SSatish Balay if (flg) { 19217699dbbSLois Curfman McInnes int rank = 0; 193d7e8b826SBarry Smith DrawLG lg; 19417699dbbSLois Curfman McInnes MPI_Initialized(&rank); 19517699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 19617699dbbSLois Curfman McInnes if (!rank) { 197*81f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 1989b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 199f8c826e1SBarry Smith snes->xmonitor = lg; 2009b94acceSBarry Smith PLogObjectParent(snes,lg); 2019b94acceSBarry Smith } 2029b94acceSBarry Smith } 2033c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2043c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2059b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2069b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 20794a424c1SBarry Smith PLogInfo(snes, 20831615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 20931615d04SLois Curfman McInnes } 21031615d04SLois Curfman McInnes else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 21131615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 21231615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 21394a424c1SBarry Smith PLogInfo(snes, 21431615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2159b94acceSBarry Smith } 216639f9d9dSBarry Smith 217639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 218639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 219639f9d9dSBarry Smith } 220639f9d9dSBarry Smith 2219b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2229b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 2239b94acceSBarry Smith if (!snes->setfromoptions) return 0; 2249b94acceSBarry Smith return (*snes->setfromoptions)(snes); 2259b94acceSBarry Smith } 2269b94acceSBarry Smith 2279b94acceSBarry Smith /*@ 2289b94acceSBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2299b94acceSBarry Smith 2309b94acceSBarry Smith Input Parameter: 2319b94acceSBarry Smith . snes - the SNES context 2329b94acceSBarry Smith 233a703fe33SLois Curfman McInnes Options Database Keys: 234a703fe33SLois Curfman McInnes $ -help, -h 235a703fe33SLois Curfman McInnes 2369b94acceSBarry Smith .keywords: SNES, nonlinear, help 2379b94acceSBarry Smith 238682d7d0cSBarry Smith .seealso: SNESSetFromOptions() 2399b94acceSBarry Smith @*/ 2409b94acceSBarry Smith int SNESPrintHelp(SNES snes) 2419b94acceSBarry Smith { 2426daaf66cSBarry Smith char p[64]; 2436daaf66cSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2446daaf66cSBarry Smith 24577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2466daaf66cSBarry Smith 2476daaf66cSBarry Smith PetscStrcpy(p,"-"); 2486daaf66cSBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2496daaf66cSBarry Smith 2506daaf66cSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2516daaf66cSBarry Smith 252d31a3109SLois Curfman McInnes PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 25333455573SSatish Balay SNESPrintTypes_Private(snes->comm,p,"snes_type"); 25477c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 255b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 256b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 257b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 258b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 259b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 260b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 2615bbad29bSBarry Smith PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 2625bbad29bSBarry Smith PetscPrintf(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 263d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor\n",p); 264d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 265b18e04deSLois Curfman McInnes if (snes->type == SNES_NONLINEAR_EQUATIONS) { 26677c4ece6SBarry Smith PetscPrintf(snes->comm, 267d31a3109SLois Curfman McInnes " Options for solving systems of nonlinear equations only:\n"); 26877c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 26977c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 270ef1dfb11SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 2711650c7b0SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p); 2721650c7b0SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p); 27377c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 27477c4ece6SBarry Smith PetscPrintf(snes->comm, 275b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 27677c4ece6SBarry Smith PetscPrintf(snes->comm, 277b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 27877c4ece6SBarry Smith PetscPrintf(snes->comm, 279b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 28077c4ece6SBarry Smith PetscPrintf(snes->comm, 281b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 28277c4ece6SBarry Smith PetscPrintf(snes->comm, 283b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 28477c4ece6SBarry Smith PetscPrintf(snes->comm, 285b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 28677c4ece6SBarry Smith PetscPrintf(snes->comm, 287b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 288b18e04deSLois Curfman McInnes } 289b18e04deSLois Curfman McInnes else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 29077c4ece6SBarry Smith PetscPrintf(snes->comm, 291d31a3109SLois Curfman McInnes " Options for solving unconstrained minimization problems only:\n"); 292b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 29377c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 29477c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 295d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p); 296d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p); 297b18e04deSLois Curfman McInnes } 2984537a946SLois Curfman McInnes PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 29977c4ece6SBarry Smith PetscPrintf(snes->comm,"a particular method\n"); 3006daaf66cSBarry Smith if (snes->printhelp) (*snes->printhelp)(snes,p); 3019b94acceSBarry Smith return 0; 3029b94acceSBarry Smith } 3039b94acceSBarry Smith /*@ 3049b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 3059b94acceSBarry Smith the nonlinear solvers. 3069b94acceSBarry Smith 3079b94acceSBarry Smith Input Parameters: 3089b94acceSBarry Smith . snes - the SNES context 3099b94acceSBarry Smith . usrP - optional user context 3109b94acceSBarry Smith 3119b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3129b94acceSBarry Smith 3139b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3149b94acceSBarry Smith @*/ 3159b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3169b94acceSBarry Smith { 31777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3189b94acceSBarry Smith snes->user = usrP; 3199b94acceSBarry Smith return 0; 3209b94acceSBarry Smith } 32174679c65SBarry Smith 3229b94acceSBarry Smith /*@C 3239b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3249b94acceSBarry Smith nonlinear solvers. 3259b94acceSBarry Smith 3269b94acceSBarry Smith Input Parameter: 3279b94acceSBarry Smith . snes - SNES context 3289b94acceSBarry Smith 3299b94acceSBarry Smith Output Parameter: 3309b94acceSBarry Smith . usrP - user context 3319b94acceSBarry Smith 3329b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3339b94acceSBarry Smith 3349b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3359b94acceSBarry Smith @*/ 3369b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3379b94acceSBarry Smith { 33877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3399b94acceSBarry Smith *usrP = snes->user; 3409b94acceSBarry Smith return 0; 3419b94acceSBarry Smith } 34274679c65SBarry Smith 3439b94acceSBarry Smith /*@ 3449b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3459b94acceSBarry Smith nonlinear solver. 3469b94acceSBarry Smith 3479b94acceSBarry Smith Input Parameter: 3489b94acceSBarry Smith . snes - SNES context 3499b94acceSBarry Smith 3509b94acceSBarry Smith Output Parameter: 3519b94acceSBarry Smith . iter - iteration number 3529b94acceSBarry Smith 3539b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3549b94acceSBarry Smith @*/ 3559b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3569b94acceSBarry Smith { 35777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 35874679c65SBarry Smith PetscValidIntPointer(iter); 3599b94acceSBarry Smith *iter = snes->iter; 3609b94acceSBarry Smith return 0; 3619b94acceSBarry Smith } 36274679c65SBarry Smith 3639b94acceSBarry Smith /*@ 3649b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3659b94acceSBarry Smith with SNESSSetFunction(). 3669b94acceSBarry Smith 3679b94acceSBarry Smith Input Parameter: 3689b94acceSBarry Smith . snes - SNES context 3699b94acceSBarry Smith 3709b94acceSBarry Smith Output Parameter: 3719b94acceSBarry Smith . fnorm - 2-norm of function 3729b94acceSBarry Smith 3739b94acceSBarry Smith Note: 3749b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 375a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 376a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 3779b94acceSBarry Smith 3789b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 379a86d99e1SLois Curfman McInnes 380a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 3819b94acceSBarry Smith @*/ 3829b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 3839b94acceSBarry Smith { 38477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 38574679c65SBarry Smith PetscValidScalarPointer(fnorm); 38674679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 38774679c65SBarry Smith SETERRQ(1,"SNESGetFunctionNorm:For SNES_NONLINEAR_EQUATIONS only"); 38874679c65SBarry Smith } 3899b94acceSBarry Smith *fnorm = snes->norm; 3909b94acceSBarry Smith return 0; 3919b94acceSBarry Smith } 39274679c65SBarry Smith 3939b94acceSBarry Smith /*@ 3949b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 3959b94acceSBarry Smith with SNESSSetGradient(). 3969b94acceSBarry Smith 3979b94acceSBarry Smith Input Parameter: 3989b94acceSBarry Smith . snes - SNES context 3999b94acceSBarry Smith 4009b94acceSBarry Smith Output Parameter: 4019b94acceSBarry Smith . fnorm - 2-norm of gradient 4029b94acceSBarry Smith 4039b94acceSBarry Smith Note: 4049b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 405a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 406a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4079b94acceSBarry Smith 4089b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 409a86d99e1SLois Curfman McInnes 410a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4119b94acceSBarry Smith @*/ 4129b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4139b94acceSBarry Smith { 41477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 41574679c65SBarry Smith PetscValidScalarPointer(gnorm); 41674679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 41774679c65SBarry Smith SETERRQ(1,"SNESGetGradientNorm:For SNES_UNCONSTRAINED_MINIMIZATION only"); 41874679c65SBarry Smith } 4199b94acceSBarry Smith *gnorm = snes->norm; 4209b94acceSBarry Smith return 0; 4219b94acceSBarry Smith } 42274679c65SBarry Smith 4239b94acceSBarry Smith /*@ 4249b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4259b94acceSBarry Smith attempted by the nonlinear solver. 4269b94acceSBarry Smith 4279b94acceSBarry Smith Input Parameter: 4289b94acceSBarry Smith . snes - SNES context 4299b94acceSBarry Smith 4309b94acceSBarry Smith Output Parameter: 4319b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4329b94acceSBarry Smith 4339b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4349b94acceSBarry Smith @*/ 4359b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4369b94acceSBarry Smith { 43777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 43874679c65SBarry Smith PetscValidIntPointer(nfails); 4399b94acceSBarry Smith *nfails = snes->nfailures; 4409b94acceSBarry Smith return 0; 4419b94acceSBarry Smith } 4429b94acceSBarry Smith /*@C 4439b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 4449b94acceSBarry Smith 4459b94acceSBarry Smith Input Parameter: 4469b94acceSBarry Smith . snes - the SNES context 4479b94acceSBarry Smith 4489b94acceSBarry Smith Output Parameter: 4499b94acceSBarry Smith . sles - the SLES context 4509b94acceSBarry Smith 4519b94acceSBarry Smith Notes: 4529b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 4539b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 4549b94acceSBarry Smith KSP and PC contexts as well. 4559b94acceSBarry Smith 4569b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 4579b94acceSBarry Smith 4589b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 4599b94acceSBarry Smith @*/ 4609b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 4619b94acceSBarry Smith { 46277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 4639b94acceSBarry Smith *sles = snes->sles; 4649b94acceSBarry Smith return 0; 4659b94acceSBarry Smith } 4669b94acceSBarry Smith /* -----------------------------------------------------------*/ 4679b94acceSBarry Smith /*@C 4689b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 4699b94acceSBarry Smith 4709b94acceSBarry Smith Input Parameter: 4719b94acceSBarry Smith . comm - MPI communicator 4729b94acceSBarry Smith . type - type of method, one of 4739b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS 4749b94acceSBarry Smith $ (for systems of nonlinear equations) 4759b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION 4769b94acceSBarry Smith $ (for unconstrained minimization) 4779b94acceSBarry Smith 4789b94acceSBarry Smith Output Parameter: 4799b94acceSBarry Smith . outsnes - the new SNES context 4809b94acceSBarry Smith 481c1f60f51SBarry Smith Options Database Key: 48219bd6747SLois Curfman McInnes $ -snes_mf - use default matrix-free Jacobian-vector products, 48319bd6747SLois Curfman McInnes $ and no preconditioning matrix 48419bd6747SLois Curfman McInnes $ -snes_mf_operator - use default matrix-free Jacobian-vector 48519bd6747SLois Curfman McInnes $ products, and a user-provided preconditioning matrix 48619bd6747SLois Curfman McInnes $ as set by SNESSetJacobian() 48719bd6747SLois Curfman McInnes $ -snes_fd - use (slow!) finite differences to compute Jacobian 488c1f60f51SBarry Smith 4899b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 4909b94acceSBarry Smith 49163a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 4929b94acceSBarry Smith @*/ 4934b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 4949b94acceSBarry Smith { 4959b94acceSBarry Smith int ierr; 4969b94acceSBarry Smith SNES snes; 4979b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 49837fcc0dbSBarry Smith 4999b94acceSBarry Smith *outsnes = 0; 5002a0acf01SLois Curfman McInnes if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 5012a0acf01SLois Curfman McInnes SETERRQ(1,"SNESCreate:incorrect method type"); 5020452661fSBarry Smith PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm); 5039b94acceSBarry Smith PLogObjectCreate(snes); 5049b94acceSBarry Smith snes->max_its = 50; 5059b94acceSBarry Smith snes->max_funcs = 1000; 5069b94acceSBarry Smith snes->norm = 0.0; 5075a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 508b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 509b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5109b94acceSBarry Smith snes->atol = 1.e-10; 5115a2d0531SBarry Smith } 512b4874afaSBarry Smith else { 513b4874afaSBarry Smith snes->rtol = 1.e-8; 514b4874afaSBarry Smith snes->ttol = 0.0; 515b4874afaSBarry Smith snes->atol = 1.e-50; 516b4874afaSBarry Smith } 5179b94acceSBarry Smith snes->xtol = 1.e-8; 518b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5199b94acceSBarry Smith snes->nfuncs = 0; 5209b94acceSBarry Smith snes->nfailures = 0; 521639f9d9dSBarry Smith snes->numbermonitors = 0; 5229b94acceSBarry Smith snes->data = 0; 5239b94acceSBarry Smith snes->view = 0; 5249b94acceSBarry Smith snes->computeumfunction = 0; 5259b94acceSBarry Smith snes->umfunP = 0; 5269b94acceSBarry Smith snes->fc = 0; 5279b94acceSBarry Smith snes->deltatol = 1.e-12; 5289b94acceSBarry Smith snes->fmin = -1.e30; 5299b94acceSBarry Smith snes->method_class = type; 5309b94acceSBarry Smith snes->set_method_called = 0; 531a703fe33SLois Curfman McInnes snes->setup_called = 0; 5329b94acceSBarry Smith snes->ksp_ewconv = 0; 5337d1a2b2bSBarry Smith snes->type = -1; 5346f24a144SLois Curfman McInnes snes->xmonitor = 0; 5356f24a144SLois Curfman McInnes snes->vwork = 0; 5366f24a144SLois Curfman McInnes snes->nwork = 0; 5379b94acceSBarry Smith 5389b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 5390452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 5409b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 5419b94acceSBarry Smith kctx->version = 2; 5429b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 5439b94acceSBarry Smith this was too large for some test cases */ 5449b94acceSBarry Smith kctx->rtol_last = 0; 5459b94acceSBarry Smith kctx->rtol_max = .9; 5469b94acceSBarry Smith kctx->gamma = 1.0; 5479b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 5489b94acceSBarry Smith kctx->alpha = kctx->alpha2; 5499b94acceSBarry Smith kctx->threshold = .1; 5509b94acceSBarry Smith kctx->lresid_last = 0; 5519b94acceSBarry Smith kctx->norm_last = 0; 5529b94acceSBarry Smith 5539b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 5549b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 5555334005bSBarry Smith 5569b94acceSBarry Smith *outsnes = snes; 5579b94acceSBarry Smith return 0; 5589b94acceSBarry Smith } 5599b94acceSBarry Smith 5609b94acceSBarry Smith /* --------------------------------------------------------------- */ 5619b94acceSBarry Smith /*@C 5629b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 5639b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 5649b94acceSBarry Smith equations. 5659b94acceSBarry Smith 5669b94acceSBarry Smith Input Parameters: 5679b94acceSBarry Smith . snes - the SNES context 5689b94acceSBarry Smith . func - function evaluation routine 5699b94acceSBarry Smith . r - vector to store function value 5702cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 5712cd2dfdcSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 5729b94acceSBarry Smith 5739b94acceSBarry Smith Calling sequence of func: 574313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Vec f,void *ctx); 5759b94acceSBarry Smith 5769b94acceSBarry Smith . x - input vector 577313e4042SLois Curfman McInnes . f - function vector 5782cd2dfdcSLois Curfman McInnes . ctx - optional user-defined function context 5799b94acceSBarry Smith 5809b94acceSBarry Smith Notes: 5819b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 5829b94acceSBarry Smith $ f'(x) x = -f(x), 5839b94acceSBarry Smith $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 5849b94acceSBarry Smith 5859b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 5869b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 5879b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 5889b94acceSBarry Smith 5899b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 5909b94acceSBarry Smith 591a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 5929b94acceSBarry Smith @*/ 5935334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 5949b94acceSBarry Smith { 59577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5969b94acceSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1, 59748d91487SBarry Smith "SNESSetFunction:For SNES_NONLINEAR_EQUATIONS only"); 5989b94acceSBarry Smith snes->computefunction = func; 5999b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6009b94acceSBarry Smith snes->funP = ctx; 6019b94acceSBarry Smith return 0; 6029b94acceSBarry Smith } 6039b94acceSBarry Smith 6049b94acceSBarry Smith /*@ 6059b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6069b94acceSBarry Smith SNESSetFunction(). 6079b94acceSBarry Smith 6089b94acceSBarry Smith Input Parameters: 6099b94acceSBarry Smith . snes - the SNES context 6109b94acceSBarry Smith . x - input vector 6119b94acceSBarry Smith 6129b94acceSBarry Smith Output Parameter: 6133638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6149b94acceSBarry Smith 6151bffabb2SLois Curfman McInnes Notes: 6169b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6179b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6189b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6199b94acceSBarry Smith 6209b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6219b94acceSBarry Smith 622a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6239b94acceSBarry Smith @*/ 6249b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 6259b94acceSBarry Smith { 6269b94acceSBarry Smith int ierr; 6279b94acceSBarry Smith 62874679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 62974679c65SBarry Smith SETERRQ(1,"SNESComputeFunction: For SNES_NONLINEAR_EQUATIONS only"); 6309b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 6319b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 6329b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 6339b94acceSBarry Smith return 0; 6349b94acceSBarry Smith } 6359b94acceSBarry Smith 6369b94acceSBarry Smith /*@C 6379b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 6389b94acceSBarry Smith unconstrained minimization. 6399b94acceSBarry Smith 6409b94acceSBarry Smith Input Parameters: 6419b94acceSBarry Smith . snes - the SNES context 6429b94acceSBarry Smith . func - function evaluation routine 643044dda88SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 644044dda88SLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6459b94acceSBarry Smith 6469b94acceSBarry Smith Calling sequence of func: 6479b94acceSBarry Smith . func (SNES snes,Vec x,double *f,void *ctx); 6489b94acceSBarry Smith 6499b94acceSBarry Smith . x - input vector 6509b94acceSBarry Smith . f - function 651044dda88SLois Curfman McInnes . ctx - [optional] user-defined function context 6529b94acceSBarry Smith 6539b94acceSBarry Smith Notes: 6549b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 6559b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 6569b94acceSBarry Smith SNESSetFunction(). 6579b94acceSBarry Smith 6589b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 6599b94acceSBarry Smith 660a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 661a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 6629b94acceSBarry Smith @*/ 6639b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 6649b94acceSBarry Smith void *ctx) 6659b94acceSBarry Smith { 66677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 6679b94acceSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 66848d91487SBarry Smith "SNESSetMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION"); 6699b94acceSBarry Smith snes->computeumfunction = func; 6709b94acceSBarry Smith snes->umfunP = ctx; 6719b94acceSBarry Smith return 0; 6729b94acceSBarry Smith } 6739b94acceSBarry Smith 6749b94acceSBarry Smith /*@ 6759b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 6769b94acceSBarry Smith set with SNESSetMinimizationFunction(). 6779b94acceSBarry Smith 6789b94acceSBarry Smith Input Parameters: 6799b94acceSBarry Smith . snes - the SNES context 6809b94acceSBarry Smith . x - input vector 6819b94acceSBarry Smith 6829b94acceSBarry Smith Output Parameter: 6839b94acceSBarry Smith . y - function value 6849b94acceSBarry Smith 6859b94acceSBarry Smith Notes: 6869b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 6879b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 6889b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 689a86d99e1SLois Curfman McInnes 690a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 691a86d99e1SLois Curfman McInnes 692a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 693a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 6949b94acceSBarry Smith @*/ 6959b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 6969b94acceSBarry Smith { 6979b94acceSBarry Smith int ierr; 6989b94acceSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 69948d91487SBarry Smith "SNESComputeMinimizationFunction:Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7009b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 7019b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 7029b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7039b94acceSBarry Smith return 0; 7049b94acceSBarry Smith } 7059b94acceSBarry Smith 7069b94acceSBarry Smith /*@C 7079b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 7089b94acceSBarry Smith vector for use by the SNES routines. 7099b94acceSBarry Smith 7109b94acceSBarry Smith Input Parameters: 7119b94acceSBarry Smith . snes - the SNES context 7129b94acceSBarry Smith . func - function evaluation routine 713044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 714044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 7159b94acceSBarry Smith . r - vector to store gradient value 7169b94acceSBarry Smith 7179b94acceSBarry Smith Calling sequence of func: 7189b94acceSBarry Smith . func (SNES, Vec x, Vec g, void *ctx); 7199b94acceSBarry Smith 7209b94acceSBarry Smith . x - input vector 7219b94acceSBarry Smith . g - gradient vector 722044dda88SLois Curfman McInnes . ctx - optional user-defined gradient context 7239b94acceSBarry Smith 7249b94acceSBarry Smith Notes: 7259b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7269b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7279b94acceSBarry Smith SNESSetFunction(). 7289b94acceSBarry Smith 7299b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 7309b94acceSBarry Smith 731a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 732a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 7339b94acceSBarry Smith @*/ 73474679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 7359b94acceSBarry Smith { 73677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 7379b94acceSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 73848d91487SBarry Smith "SNESSetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 7399b94acceSBarry Smith snes->computefunction = func; 7409b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7419b94acceSBarry Smith snes->funP = ctx; 7429b94acceSBarry Smith return 0; 7439b94acceSBarry Smith } 7449b94acceSBarry Smith 7459b94acceSBarry Smith /*@ 746a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 747a86d99e1SLois Curfman McInnes SNESSetGradient(). 7489b94acceSBarry Smith 7499b94acceSBarry Smith Input Parameters: 7509b94acceSBarry Smith . snes - the SNES context 7519b94acceSBarry Smith . x - input vector 7529b94acceSBarry Smith 7539b94acceSBarry Smith Output Parameter: 7549b94acceSBarry Smith . y - gradient vector 7559b94acceSBarry Smith 7569b94acceSBarry Smith Notes: 7579b94acceSBarry Smith SNESComputeGradient() is valid only for 7589b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7599b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 760a86d99e1SLois Curfman McInnes 761a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 762a86d99e1SLois Curfman McInnes 763a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 764a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 7659b94acceSBarry Smith @*/ 7669b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 7679b94acceSBarry Smith { 7689b94acceSBarry Smith int ierr; 76974679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 77074679c65SBarry Smith SETERRQ(1,"SNESComputeGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 7719b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 7729b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 7739b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 7749b94acceSBarry Smith return 0; 7759b94acceSBarry Smith } 7769b94acceSBarry Smith 77762fef451SLois Curfman McInnes /*@ 77862fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 77962fef451SLois Curfman McInnes set with SNESSetJacobian(). 78062fef451SLois Curfman McInnes 78162fef451SLois Curfman McInnes Input Parameters: 78262fef451SLois Curfman McInnes . snes - the SNES context 78362fef451SLois Curfman McInnes . x - input vector 78462fef451SLois Curfman McInnes 78562fef451SLois Curfman McInnes Output Parameters: 78662fef451SLois Curfman McInnes . A - Jacobian matrix 78762fef451SLois Curfman McInnes . B - optional preconditioning matrix 78862fef451SLois Curfman McInnes . flag - flag indicating matrix structure 78962fef451SLois Curfman McInnes 79062fef451SLois Curfman McInnes Notes: 79162fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 79262fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 79362fef451SLois Curfman McInnes 794dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 795dc5a77f8SLois Curfman McInnes flag parameter. 79662fef451SLois Curfman McInnes 79762fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 79862fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 79962fef451SLois Curfman McInnes methods is SNESComputeJacobian(). 80062fef451SLois Curfman McInnes 80162fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 80262fef451SLois Curfman McInnes 80362fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 80462fef451SLois Curfman McInnes @*/ 8059b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8069b94acceSBarry Smith { 8079b94acceSBarry Smith int ierr; 80874679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 80974679c65SBarry Smith SETERRQ(1,"SNESComputeJacobian: For SNES_NONLINEAR_EQUATIONS only"); 8109b94acceSBarry Smith if (!snes->computejacobian) return 0; 8119b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 812c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 8139b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 8149b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 8156d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 81677c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 81777c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 8189b94acceSBarry Smith return 0; 8199b94acceSBarry Smith } 8209b94acceSBarry Smith 82162fef451SLois Curfman McInnes /*@ 82262fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 82362fef451SLois Curfman McInnes set with SNESSetHessian(). 82462fef451SLois Curfman McInnes 82562fef451SLois Curfman McInnes Input Parameters: 82662fef451SLois Curfman McInnes . snes - the SNES context 82762fef451SLois Curfman McInnes . x - input vector 82862fef451SLois Curfman McInnes 82962fef451SLois Curfman McInnes Output Parameters: 83062fef451SLois Curfman McInnes . A - Hessian matrix 83162fef451SLois Curfman McInnes . B - optional preconditioning matrix 83262fef451SLois Curfman McInnes . flag - flag indicating matrix structure 83362fef451SLois Curfman McInnes 83462fef451SLois Curfman McInnes Notes: 83562fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 83662fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 83762fef451SLois Curfman McInnes 838dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 839dc5a77f8SLois Curfman McInnes flag parameter. 84062fef451SLois Curfman McInnes 84162fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 84262fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 84362fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 84462fef451SLois Curfman McInnes 84562fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 84662fef451SLois Curfman McInnes 847a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 848a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 84962fef451SLois Curfman McInnes @*/ 85062fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 8519b94acceSBarry Smith { 8529b94acceSBarry Smith int ierr; 85374679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 85474679c65SBarry Smith SETERRQ(1,"SNESComputeHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 8559b94acceSBarry Smith if (!snes->computejacobian) return 0; 85662fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 8570f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 85862fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 85962fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 8600f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 86177c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 86277c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 8639b94acceSBarry Smith return 0; 8649b94acceSBarry Smith } 8659b94acceSBarry Smith 8669b94acceSBarry Smith /*@C 8679b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 868044dda88SLois Curfman McInnes location to store the matrix. 8699b94acceSBarry Smith 8709b94acceSBarry Smith Input Parameters: 8719b94acceSBarry Smith . snes - the SNES context 8729b94acceSBarry Smith . A - Jacobian matrix 8739b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 8749b94acceSBarry Smith . func - Jacobian evaluation routine 8752cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 8762cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 8779b94acceSBarry Smith 8789b94acceSBarry Smith Calling sequence of func: 879313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 8809b94acceSBarry Smith 8819b94acceSBarry Smith . x - input vector 8829b94acceSBarry Smith . A - Jacobian matrix 8839b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 884ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 885ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 8862cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Jacobian context 8879b94acceSBarry Smith 8889b94acceSBarry Smith Notes: 889dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 8902cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 891ac21db08SLois Curfman McInnes 892ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 8939b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 8949b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 8959b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 8969b94acceSBarry Smith throughout the global iterations. 8979b94acceSBarry Smith 8989b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 8999b94acceSBarry Smith 900ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 9019b94acceSBarry Smith @*/ 9029b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 9039b94acceSBarry Smith MatStructure*,void*),void *ctx) 9049b94acceSBarry Smith { 90577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 90674679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 90774679c65SBarry Smith SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only"); 9089b94acceSBarry Smith snes->computejacobian = func; 9099b94acceSBarry Smith snes->jacP = ctx; 9109b94acceSBarry Smith snes->jacobian = A; 9119b94acceSBarry Smith snes->jacobian_pre = B; 9129b94acceSBarry Smith return 0; 9139b94acceSBarry Smith } 91462fef451SLois Curfman McInnes 915b4fd4287SBarry Smith /*@ 916b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 917b4fd4287SBarry Smith provided context for evaluating the Jacobian. 918b4fd4287SBarry Smith 919b4fd4287SBarry Smith Input Parameter: 920b4fd4287SBarry Smith . snes - the nonlinear solver context 921b4fd4287SBarry Smith 922b4fd4287SBarry Smith Output Parameters: 923b4fd4287SBarry Smith . A - location to stash Jacobian matrix (or PETSC_NULL) 924b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 925b4fd4287SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 926b4fd4287SBarry Smith 927b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 928b4fd4287SBarry Smith @*/ 929b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 930b4fd4287SBarry Smith { 93174679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 93274679c65SBarry Smith SETERRQ(1,"SNESSetJacobian:For SNES_NONLINEAR_EQUATIONS only"); 933b4fd4287SBarry Smith if (A) *A = snes->jacobian; 934b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 935b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 936b4fd4287SBarry Smith return 0; 937b4fd4287SBarry Smith } 938b4fd4287SBarry Smith 9399b94acceSBarry Smith /*@C 9409b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 941044dda88SLois Curfman McInnes location to store the matrix. 9429b94acceSBarry Smith 9439b94acceSBarry Smith Input Parameters: 9449b94acceSBarry Smith . snes - the SNES context 9459b94acceSBarry Smith . A - Hessian matrix 9469b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 9479b94acceSBarry Smith . func - Jacobian evaluation routine 948313e4042SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 949313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 9509b94acceSBarry Smith 9519b94acceSBarry Smith Calling sequence of func: 952313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9539b94acceSBarry Smith 9549b94acceSBarry Smith . x - input vector 9559b94acceSBarry Smith . A - Hessian matrix 9569b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 957ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 958ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 9592cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Hessian context 9609b94acceSBarry Smith 9619b94acceSBarry Smith Notes: 962dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 9632cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 964ac21db08SLois Curfman McInnes 9659b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 9669b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 9679b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9689b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9699b94acceSBarry Smith throughout the global iterations. 9709b94acceSBarry Smith 9719b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 9729b94acceSBarry Smith 973ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 9749b94acceSBarry Smith @*/ 9759b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 9769b94acceSBarry Smith MatStructure*,void*),void *ctx) 9779b94acceSBarry Smith { 97877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 97974679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 98074679c65SBarry Smith SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 9819b94acceSBarry Smith snes->computejacobian = func; 9829b94acceSBarry Smith snes->jacP = ctx; 9839b94acceSBarry Smith snes->jacobian = A; 9849b94acceSBarry Smith snes->jacobian_pre = B; 9859b94acceSBarry Smith return 0; 9869b94acceSBarry Smith } 9879b94acceSBarry Smith 98862fef451SLois Curfman McInnes /*@ 98962fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 99062fef451SLois Curfman McInnes provided context for evaluating the Hessian. 99162fef451SLois Curfman McInnes 99262fef451SLois Curfman McInnes Input Parameter: 99362fef451SLois Curfman McInnes . snes - the nonlinear solver context 99462fef451SLois Curfman McInnes 99562fef451SLois Curfman McInnes Output Parameters: 99662fef451SLois Curfman McInnes . A - location to stash Hessian matrix (or PETSC_NULL) 99762fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 99862fef451SLois Curfman McInnes . ctx - location to stash Hessian ctx (or PETSC_NULL) 99962fef451SLois Curfman McInnes 100062fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 100162fef451SLois Curfman McInnes @*/ 100262fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 100362fef451SLois Curfman McInnes { 100474679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 100574679c65SBarry Smith SETERRQ(1,"SNESSetHessian:For SNES_UNCONSTRAINED_MINIMIZATION only"); 100662fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 100762fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 100862fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 100962fef451SLois Curfman McInnes return 0; 101062fef451SLois Curfman McInnes } 101162fef451SLois Curfman McInnes 10129b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 10139b94acceSBarry Smith 10149b94acceSBarry Smith /*@ 10159b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1016272ac6f2SLois Curfman McInnes of a nonlinear solver. 10179b94acceSBarry Smith 10189b94acceSBarry Smith Input Parameter: 10199b94acceSBarry Smith . snes - the SNES context 10208ddd3da0SLois Curfman McInnes . x - the solution vector 10219b94acceSBarry Smith 1022272ac6f2SLois Curfman McInnes Notes: 1023272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1024272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1025272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1026272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1027272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1028272ac6f2SLois Curfman McInnes 10299b94acceSBarry Smith .keywords: SNES, nonlinear, setup 10309b94acceSBarry Smith 10319b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 10329b94acceSBarry Smith @*/ 10338ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 10349b94acceSBarry Smith { 1035272ac6f2SLois Curfman McInnes int ierr, flg; 103677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 103777c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 10388ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 10399b94acceSBarry Smith 1040c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1041c1f60f51SBarry Smith /* 1042c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1043dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1044c1f60f51SBarry Smith */ 1045c1f60f51SBarry Smith if (flg) { 1046c1f60f51SBarry Smith Mat J; 1047c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1048c1f60f51SBarry Smith PLogObjectParent(snes,J); 1049c1f60f51SBarry Smith snes->mfshell = J; 1050c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1051c1f60f51SBarry Smith snes->jacobian = J; 1052c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1053c1f60f51SBarry Smith } 1054c1f60f51SBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1055c1f60f51SBarry Smith snes->jacobian = J; 1056c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1057c1f60f51SBarry Smith } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free operator option"); 1058c1f60f51SBarry Smith } 1059272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1060c1f60f51SBarry Smith /* 1061dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1062c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1063c1f60f51SBarry Smith */ 106431615d04SLois Curfman McInnes if (flg) { 1065272ac6f2SLois Curfman McInnes Mat J; 1066272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1067272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1068272ac6f2SLois Curfman McInnes snes->mfshell = J; 106983e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 107083e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 107194a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 107283e56afdSLois Curfman McInnes } 107383e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 107483e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 107594a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 107683e56afdSLois Curfman McInnes } else SETERRQ(1,"SNESSetUp:Method class doesn't support matrix-free option"); 1077272ac6f2SLois Curfman McInnes } 10789b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 107937fcc0dbSBarry Smith if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first"); 108037fcc0dbSBarry Smith if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetFunction() first"); 108137fcc0dbSBarry Smith if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetJacobian() first"); 1082d804570eSBarry Smith if (snes->vec_func == snes->vec_sol) SETERRQ(1,"SNESSetUp:Solution vector cannot be function vector"); 1083a703fe33SLois Curfman McInnes 1084a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 108540191667SLois Curfman McInnes if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1086a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1087a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1088a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1089a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1090a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1091a703fe33SLois Curfman McInnes } 10929b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 109337fcc0dbSBarry Smith if (!snes->vec_func) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first"); 109437fcc0dbSBarry Smith if (!snes->computefunction) SETERRQ(1,"SNESSetUp:Must call SNESSetGradient() first"); 109537fcc0dbSBarry Smith if (!snes->computeumfunction) 109637fcc0dbSBarry Smith SETERRQ(1,"SNESSetUp:Must call SNESSetMinimizationFunction() first"); 109737fcc0dbSBarry Smith if (!snes->jacobian) SETERRQ(1,"SNESSetUp:Must call SNESSetHessian() first"); 10989b94acceSBarry Smith } else SETERRQ(1,"SNESSetUp:Unknown method class"); 1099a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1100a703fe33SLois Curfman McInnes snes->setup_called = 1; 1101a703fe33SLois Curfman McInnes return 0; 11029b94acceSBarry Smith } 11039b94acceSBarry Smith 11049b94acceSBarry Smith /*@C 11059b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 11069b94acceSBarry Smith with SNESCreate(). 11079b94acceSBarry Smith 11089b94acceSBarry Smith Input Parameter: 11099b94acceSBarry Smith . snes - the SNES context 11109b94acceSBarry Smith 11119b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 11129b94acceSBarry Smith 111363a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 11149b94acceSBarry Smith @*/ 11159b94acceSBarry Smith int SNESDestroy(SNES snes) 11169b94acceSBarry Smith { 11179b94acceSBarry Smith int ierr; 111877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 11199b94acceSBarry Smith ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr); 11200452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 11219b94acceSBarry Smith if (snes->mfshell) MatDestroy(snes->mfshell); 11229b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 11233e2e452bSBarry Smith if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 11246f24a144SLois Curfman McInnes if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 11259b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 11260452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 11279b94acceSBarry Smith return 0; 11289b94acceSBarry Smith } 11299b94acceSBarry Smith 11309b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 11319b94acceSBarry Smith 11329b94acceSBarry Smith /*@ 1133d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 11349b94acceSBarry Smith 11359b94acceSBarry Smith Input Parameters: 11369b94acceSBarry Smith . snes - the SNES context 113733174efeSLois Curfman McInnes . atol - absolute convergence tolerance 113833174efeSLois Curfman McInnes . rtol - relative convergence tolerance 113933174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 114033174efeSLois Curfman McInnes of the change in the solution between steps 114133174efeSLois Curfman McInnes . maxit - maximum number of iterations 114233174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 11439b94acceSBarry Smith 114433174efeSLois Curfman McInnes Options Database Keys: 114533174efeSLois Curfman McInnes $ -snes_atol <atol> 114633174efeSLois Curfman McInnes $ -snes_rtol <rtol> 114733174efeSLois Curfman McInnes $ -snes_stol <stol> 114833174efeSLois Curfman McInnes $ -snes_max_it <maxit> 114933174efeSLois Curfman McInnes $ -snes_max_funcs <maxf> 11509b94acceSBarry Smith 1151d7a720efSLois Curfman McInnes Notes: 11529b94acceSBarry Smith The default maximum number of iterations is 50. 11539b94acceSBarry Smith The default maximum number of function evaluations is 1000. 11549b94acceSBarry Smith 115533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 11569b94acceSBarry Smith 1157d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 11589b94acceSBarry Smith @*/ 1159d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 11609b94acceSBarry Smith { 116177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1162d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1163d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1164d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1165d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1166d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 11679b94acceSBarry Smith return 0; 11689b94acceSBarry Smith } 11699b94acceSBarry Smith 11709b94acceSBarry Smith /*@ 117133174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 117233174efeSLois Curfman McInnes 117333174efeSLois Curfman McInnes Input Parameters: 117433174efeSLois Curfman McInnes . snes - the SNES context 117533174efeSLois Curfman McInnes . atol - absolute convergence tolerance 117633174efeSLois Curfman McInnes . rtol - relative convergence tolerance 117733174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 117833174efeSLois Curfman McInnes of the change in the solution between steps 117933174efeSLois Curfman McInnes . maxit - maximum number of iterations 118033174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 118133174efeSLois Curfman McInnes 118233174efeSLois Curfman McInnes Notes: 118333174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 118433174efeSLois Curfman McInnes 118533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 118633174efeSLois Curfman McInnes 118733174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 118833174efeSLois Curfman McInnes @*/ 118933174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 119033174efeSLois Curfman McInnes { 119133174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 119233174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 119333174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 119433174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 119533174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 119633174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 119733174efeSLois Curfman McInnes return 0; 119833174efeSLois Curfman McInnes } 119933174efeSLois Curfman McInnes 120033174efeSLois Curfman McInnes /*@ 12019b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 12029b94acceSBarry Smith 12039b94acceSBarry Smith Input Parameters: 12049b94acceSBarry Smith . snes - the SNES context 12059b94acceSBarry Smith . tol - tolerance 12069b94acceSBarry Smith 12079b94acceSBarry Smith Options Database Key: 1208d7a720efSLois Curfman McInnes $ -snes_trtol <tol> 12099b94acceSBarry Smith 12109b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 12119b94acceSBarry Smith 1212d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 12139b94acceSBarry Smith @*/ 12149b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 12159b94acceSBarry Smith { 121677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12179b94acceSBarry Smith snes->deltatol = tol; 12189b94acceSBarry Smith return 0; 12199b94acceSBarry Smith } 12209b94acceSBarry Smith 12219b94acceSBarry Smith /*@ 122277c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 12239b94acceSBarry Smith for unconstrained minimization solvers. 12249b94acceSBarry Smith 12259b94acceSBarry Smith Input Parameters: 12269b94acceSBarry Smith . snes - the SNES context 12279b94acceSBarry Smith . ftol - minimum function tolerance 12289b94acceSBarry Smith 12299b94acceSBarry Smith Options Database Key: 1230d7a720efSLois Curfman McInnes $ -snes_fmin <ftol> 12319b94acceSBarry Smith 12329b94acceSBarry Smith Note: 123377c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 12349b94acceSBarry Smith methods only. 12359b94acceSBarry Smith 12369b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 12379b94acceSBarry Smith 1238d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 12399b94acceSBarry Smith @*/ 124077c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 12419b94acceSBarry Smith { 124277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12439b94acceSBarry Smith snes->fmin = ftol; 12449b94acceSBarry Smith return 0; 12459b94acceSBarry Smith } 12469b94acceSBarry Smith 12479b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 12489b94acceSBarry Smith 12499b94acceSBarry Smith /*@C 12509b94acceSBarry Smith SNESSetMonitor - Sets the function that is to be used at every 12519b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 12529b94acceSBarry Smith progress. 12539b94acceSBarry Smith 12549b94acceSBarry Smith Input Parameters: 12559b94acceSBarry Smith . snes - the SNES context 12569b94acceSBarry Smith . func - monitoring routine 1257044dda88SLois Curfman McInnes . mctx - [optional] user-defined context for private data for the 1258044dda88SLois Curfman McInnes monitor routine (may be PETSC_NULL) 12599b94acceSBarry Smith 12609b94acceSBarry Smith Calling sequence of func: 1261682d7d0cSBarry Smith int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 12629b94acceSBarry Smith 12639b94acceSBarry Smith $ snes - the SNES context 12649b94acceSBarry Smith $ its - iteration number 1265044dda88SLois Curfman McInnes $ mctx - [optional] monitoring context 12669b94acceSBarry Smith $ 12679b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 12689b94acceSBarry Smith $ norm - 2-norm function value (may be estimated) 12699b94acceSBarry Smith $ 12709b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 12719b94acceSBarry Smith $ norm - 2-norm gradient value (may be estimated) 12729b94acceSBarry Smith 12739665c990SLois Curfman McInnes Options Database Keys: 12749665c990SLois Curfman McInnes $ -snes_monitor : sets SNESDefaultMonitor() 12759665c990SLois Curfman McInnes $ -snes_xmonitor : sets line graph monitor, 12769665c990SLois Curfman McInnes $ uses SNESLGMonitorCreate() 12779665c990SLois Curfman McInnes $ -snes_cancelmonitors : cancels all monitors that have 12789665c990SLois Curfman McInnes $ been hardwired into a code by 12799665c990SLois Curfman McInnes $ calls to SNESSetMonitor(), but 12809665c990SLois Curfman McInnes $ does not cancel those set via 12819665c990SLois Curfman McInnes $ the options database. 12829665c990SLois Curfman McInnes 12839665c990SLois Curfman McInnes 1284639f9d9dSBarry Smith Notes: 12856bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 12866bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 12876bc08f3fSLois Curfman McInnes order in which they were set. 1288639f9d9dSBarry Smith 12899b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 12909b94acceSBarry Smith 12919b94acceSBarry Smith .seealso: SNESDefaultMonitor() 12929b94acceSBarry Smith @*/ 129374679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 12949b94acceSBarry Smith { 1295639f9d9dSBarry Smith if (!func) { 1296639f9d9dSBarry Smith snes->numbermonitors = 0; 1297639f9d9dSBarry Smith return 0; 1298639f9d9dSBarry Smith } 1299639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1300639f9d9dSBarry Smith SETERRQ(1,"SNESSetMonitor:Too many monitors set"); 1301639f9d9dSBarry Smith } 1302639f9d9dSBarry Smith 1303639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1304639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 13059b94acceSBarry Smith return 0; 13069b94acceSBarry Smith } 13079b94acceSBarry Smith 13089b94acceSBarry Smith /*@C 13099b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13109b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13119b94acceSBarry Smith 13129b94acceSBarry Smith Input Parameters: 13139b94acceSBarry Smith . snes - the SNES context 13149b94acceSBarry Smith . func - routine to test for convergence 1315044dda88SLois Curfman McInnes . cctx - [optional] context for private data for the convergence routine 1316044dda88SLois Curfman McInnes (may be PETSC_NULL) 13179b94acceSBarry Smith 13189b94acceSBarry Smith Calling sequence of func: 13199b94acceSBarry Smith int func (SNES snes,double xnorm,double gnorm, 13209b94acceSBarry Smith double f,void *cctx) 13219b94acceSBarry Smith 13229b94acceSBarry Smith $ snes - the SNES context 1323044dda88SLois Curfman McInnes $ cctx - [optional] convergence context 13249b94acceSBarry Smith $ xnorm - 2-norm of current iterate 13259b94acceSBarry Smith $ 13269b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13279b94acceSBarry Smith $ gnorm - 2-norm of current step 13289b94acceSBarry Smith $ f - 2-norm of function 13299b94acceSBarry Smith $ 13309b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13319b94acceSBarry Smith $ gnorm - 2-norm of current gradient 13329b94acceSBarry Smith $ f - function value 13339b94acceSBarry Smith 13349b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13359b94acceSBarry Smith 133640191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 133740191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 13389b94acceSBarry Smith @*/ 133974679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 13409b94acceSBarry Smith { 13419b94acceSBarry Smith (snes)->converged = func; 13429b94acceSBarry Smith (snes)->cnvP = cctx; 13439b94acceSBarry Smith return 0; 13449b94acceSBarry Smith } 13459b94acceSBarry Smith 13469b94acceSBarry Smith /* 13479b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 13489b94acceSBarry Smith positive parameter delta. 13499b94acceSBarry Smith 13509b94acceSBarry Smith Input Parameters: 13519b94acceSBarry Smith . snes - the SNES context 13529b94acceSBarry Smith . y - approximate solution of linear system 13539b94acceSBarry Smith . fnorm - 2-norm of current function 13549b94acceSBarry Smith . delta - trust region size 13559b94acceSBarry Smith 13569b94acceSBarry Smith Output Parameters: 13579b94acceSBarry Smith . gpnorm - predicted function norm at the new point, assuming local 13589b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 13599b94acceSBarry Smith region, and exceeds zero otherwise. 13609b94acceSBarry Smith . ynorm - 2-norm of the step 13619b94acceSBarry Smith 13629b94acceSBarry Smith Note: 136340191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 13649b94acceSBarry Smith is set to be the maximum allowable step size. 13659b94acceSBarry Smith 13669b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 13679b94acceSBarry Smith */ 13689b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 13699b94acceSBarry Smith double *gpnorm,double *ynorm) 13709b94acceSBarry Smith { 13719b94acceSBarry Smith double norm; 13729b94acceSBarry Smith Scalar cnorm; 1373cddf8d76SBarry Smith VecNorm(y,NORM_2, &norm ); 13749b94acceSBarry Smith if (norm > *delta) { 13759b94acceSBarry Smith norm = *delta/norm; 13769b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 13779b94acceSBarry Smith cnorm = norm; 13789b94acceSBarry Smith VecScale( &cnorm, y ); 13799b94acceSBarry Smith *ynorm = *delta; 13809b94acceSBarry Smith } else { 13819b94acceSBarry Smith *gpnorm = 0.0; 13829b94acceSBarry Smith *ynorm = norm; 13839b94acceSBarry Smith } 13849b94acceSBarry Smith return 0; 13859b94acceSBarry Smith } 13869b94acceSBarry Smith 13879b94acceSBarry Smith /*@ 13889b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 138963a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 13909b94acceSBarry Smith 13919b94acceSBarry Smith Input Parameter: 13929b94acceSBarry Smith . snes - the SNES context 13938ddd3da0SLois Curfman McInnes . x - the solution vector 13949b94acceSBarry Smith 13959b94acceSBarry Smith Output Parameter: 13969b94acceSBarry Smith its - number of iterations until termination 13979b94acceSBarry Smith 13988ddd3da0SLois Curfman McInnes Note: 13998ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 14008ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 14018ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 14028ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 14038ddd3da0SLois Curfman McInnes 14049b94acceSBarry Smith .keywords: SNES, nonlinear, solve 14059b94acceSBarry Smith 140663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 14079b94acceSBarry Smith @*/ 14088ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 14099b94acceSBarry Smith { 14103c7409f5SSatish Balay int ierr, flg; 1411052efed2SBarry Smith 141277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 141374679c65SBarry Smith PetscValidIntPointer(its); 1414c4fc05e7SBarry Smith if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1415c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 14169b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 14179b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 14189b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 14193c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 14206d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 14219b94acceSBarry Smith return 0; 14229b94acceSBarry Smith } 14239b94acceSBarry Smith 14249b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 142537fcc0dbSBarry Smith static NRList *__SNESList = 0; 14269b94acceSBarry Smith 14279b94acceSBarry Smith /*@ 14284b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 14299b94acceSBarry Smith 14309b94acceSBarry Smith Input Parameters: 14319b94acceSBarry Smith . snes - the SNES context 14329b94acceSBarry Smith . method - a known method 14339b94acceSBarry Smith 1434ae12b187SLois Curfman McInnes Options Database Command: 1435ae12b187SLois Curfman McInnes $ -snes_type <method> 1436ae12b187SLois Curfman McInnes $ Use -help for a list of available methods 1437ae12b187SLois Curfman McInnes $ (for instance, ls or tr) 1438ae12b187SLois Curfman McInnes 14399b94acceSBarry Smith Notes: 14409b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 14419b94acceSBarry Smith $ Systems of nonlinear equations: 144240191667SLois Curfman McInnes $ SNES_EQ_LS - Newton's method with line search 144340191667SLois Curfman McInnes $ SNES_EQ_TR - Newton's method with trust region 14449b94acceSBarry Smith $ Unconstrained minimization: 144540191667SLois Curfman McInnes $ SNES_UM_TR - Newton's method with trust region 144640191667SLois Curfman McInnes $ SNES_UM_LS - Newton's method with line search 14479b94acceSBarry Smith 1448ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1449ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1450ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1451ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1452ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1453ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1454ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1455ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1456ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1457ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1458a703fe33SLois Curfman McInnes 1459f536c699SSatish Balay .keywords: SNES, set, method 14609b94acceSBarry Smith @*/ 14614b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 14629b94acceSBarry Smith { 14639b94acceSBarry Smith int (*r)(SNES); 146477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 14657d1a2b2bSBarry Smith if (snes->type == (int) method) return 0; 14667d1a2b2bSBarry Smith 14679b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 146837fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 146937fcc0dbSBarry Smith if (!__SNESList) {SETERRQ(1,"SNESSetType:Could not get methods");} 147037fcc0dbSBarry Smith r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 14714b0e389bSBarry Smith if (!r) {SETERRQ(1,"SNESSetType:Unknown method");} 14720452661fSBarry Smith if (snes->data) PetscFree(snes->data); 14739b94acceSBarry Smith snes->set_method_called = 1; 14749b94acceSBarry Smith return (*r)(snes); 14759b94acceSBarry Smith } 14769b94acceSBarry Smith 14779b94acceSBarry Smith /* --------------------------------------------------------------------- */ 14789b94acceSBarry Smith /*@C 14799b94acceSBarry Smith SNESRegister - Adds the method to the nonlinear solver package, given 14804b0e389bSBarry Smith a function pointer and a nonlinear solver name of the type SNESType. 14819b94acceSBarry Smith 14829b94acceSBarry Smith Input Parameters: 148340191667SLois Curfman McInnes . name - for instance SNES_EQ_LS, SNES_EQ_TR, ... 148440191667SLois Curfman McInnes . sname - corresponding string for name 14859b94acceSBarry Smith . create - routine to create method context 14869b94acceSBarry Smith 14879b94acceSBarry Smith .keywords: SNES, nonlinear, register 14889b94acceSBarry Smith 14899b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy() 14909b94acceSBarry Smith @*/ 14919b94acceSBarry Smith int SNESRegister(int name, char *sname, int (*create)(SNES)) 14929b94acceSBarry Smith { 14939b94acceSBarry Smith int ierr; 149437fcc0dbSBarry Smith if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 149537fcc0dbSBarry Smith NRRegister( __SNESList, name, sname, (int (*)(void*))create ); 14969b94acceSBarry Smith return 0; 14979b94acceSBarry Smith } 14989b94acceSBarry Smith /* --------------------------------------------------------------------- */ 14999b94acceSBarry Smith /*@C 15009b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 15019b94acceSBarry Smith registered by SNESRegister(). 15029b94acceSBarry Smith 15039b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 15049b94acceSBarry Smith 15059b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 15069b94acceSBarry Smith @*/ 15079b94acceSBarry Smith int SNESRegisterDestroy() 15089b94acceSBarry Smith { 150937fcc0dbSBarry Smith if (__SNESList) { 151037fcc0dbSBarry Smith NRDestroy( __SNESList ); 151137fcc0dbSBarry Smith __SNESList = 0; 15129b94acceSBarry Smith } 15139b94acceSBarry Smith return 0; 15149b94acceSBarry Smith } 15159b94acceSBarry Smith 15169b94acceSBarry Smith /* 15174b0e389bSBarry Smith SNESGetTypeFromOptions_Private - Sets the selected method from the 15189b94acceSBarry Smith options database. 15199b94acceSBarry Smith 15209b94acceSBarry Smith Input Parameter: 15219b94acceSBarry Smith . ctx - the SNES context 15229b94acceSBarry Smith 15239b94acceSBarry Smith Output Parameter: 15249b94acceSBarry Smith . method - solver method 15259b94acceSBarry Smith 15269b94acceSBarry Smith Returns: 15279b94acceSBarry Smith Returns 1 if the method is found; 0 otherwise. 15289b94acceSBarry Smith 15299b94acceSBarry Smith Options Database Key: 15304b0e389bSBarry Smith $ -snes_type method 15319b94acceSBarry Smith */ 1532052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 15339b94acceSBarry Smith { 1534052efed2SBarry Smith int ierr; 15359b94acceSBarry Smith char sbuf[50]; 15365baf8537SBarry Smith 1537052efed2SBarry Smith ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1538052efed2SBarry Smith if (*flg) { 1539052efed2SBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 154037fcc0dbSBarry Smith *method = (SNESType)NRFindID( __SNESList, sbuf ); 15419b94acceSBarry Smith } 15429b94acceSBarry Smith return 0; 15439b94acceSBarry Smith } 15449b94acceSBarry Smith 15459b94acceSBarry Smith /*@C 15469a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 15479b94acceSBarry Smith 15489b94acceSBarry Smith Input Parameter: 15494b0e389bSBarry Smith . snes - nonlinear solver context 15509b94acceSBarry Smith 15519b94acceSBarry Smith Output Parameter: 15529a28b0a6SLois Curfman McInnes . method - SNES method (or use PETSC_NULL) 15539a28b0a6SLois Curfman McInnes . name - name of SNES method (or use PETSC_NULL) 15549b94acceSBarry Smith 15559b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 15569b94acceSBarry Smith @*/ 15574b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name) 15589b94acceSBarry Smith { 155906a9b0adSLois Curfman McInnes int ierr; 156037fcc0dbSBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 15614b0e389bSBarry Smith if (method) *method = (SNESType) snes->type; 156237fcc0dbSBarry Smith if (name) *name = NRFindName( __SNESList, (int) snes->type ); 15639b94acceSBarry Smith return 0; 15649b94acceSBarry Smith } 15659b94acceSBarry Smith 15669b94acceSBarry Smith #include <stdio.h> 15679b94acceSBarry Smith /* 15684b0e389bSBarry Smith SNESPrintTypes_Private - Prints the SNES methods available from the 15699b94acceSBarry Smith options database. 15709b94acceSBarry Smith 15719b94acceSBarry Smith Input Parameters: 157233455573SSatish Balay . comm - communicator (usually MPI_COMM_WORLD) 15739b94acceSBarry Smith . prefix - prefix (usually "-") 15744b0e389bSBarry Smith . name - the options database name (by default "snes_type") 15759b94acceSBarry Smith */ 157633455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 15779b94acceSBarry Smith { 15789b94acceSBarry Smith FuncList *entry; 157937fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 158037fcc0dbSBarry Smith entry = __SNESList->head; 158177c4ece6SBarry Smith PetscPrintf(comm," %s%s (one of)",prefix,name); 15829b94acceSBarry Smith while (entry) { 158377c4ece6SBarry Smith PetscPrintf(comm," %s",entry->name); 15849b94acceSBarry Smith entry = entry->next; 15859b94acceSBarry Smith } 158677c4ece6SBarry Smith PetscPrintf(comm,"\n"); 15879b94acceSBarry Smith return 0; 15889b94acceSBarry Smith } 15899b94acceSBarry Smith 15909b94acceSBarry Smith /*@C 15919b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 15929b94acceSBarry Smith stored. 15939b94acceSBarry Smith 15949b94acceSBarry Smith Input Parameter: 15959b94acceSBarry Smith . snes - the SNES context 15969b94acceSBarry Smith 15979b94acceSBarry Smith Output Parameter: 15989b94acceSBarry Smith . x - the solution 15999b94acceSBarry Smith 16009b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 16019b94acceSBarry Smith 1602a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 16039b94acceSBarry Smith @*/ 16049b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 16059b94acceSBarry Smith { 160677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16079b94acceSBarry Smith *x = snes->vec_sol_always; 16089b94acceSBarry Smith return 0; 16099b94acceSBarry Smith } 16109b94acceSBarry Smith 16119b94acceSBarry Smith /*@C 16129b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 16139b94acceSBarry Smith stored. 16149b94acceSBarry Smith 16159b94acceSBarry Smith Input Parameter: 16169b94acceSBarry Smith . snes - the SNES context 16179b94acceSBarry Smith 16189b94acceSBarry Smith Output Parameter: 16199b94acceSBarry Smith . x - the solution update 16209b94acceSBarry Smith 16219b94acceSBarry Smith Notes: 16229b94acceSBarry Smith This vector is implementation dependent. 16239b94acceSBarry Smith 16249b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 16259b94acceSBarry Smith 16269b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 16279b94acceSBarry Smith @*/ 16289b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 16299b94acceSBarry Smith { 163077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16319b94acceSBarry Smith *x = snes->vec_sol_update_always; 16329b94acceSBarry Smith return 0; 16339b94acceSBarry Smith } 16349b94acceSBarry Smith 16359b94acceSBarry Smith /*@C 16363638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 16379b94acceSBarry Smith 16389b94acceSBarry Smith Input Parameter: 16399b94acceSBarry Smith . snes - the SNES context 16409b94acceSBarry Smith 16419b94acceSBarry Smith Output Parameter: 16423638b69dSLois Curfman McInnes . r - the function 16439b94acceSBarry Smith 16449b94acceSBarry Smith Notes: 16459b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 16469b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 16479b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 16489b94acceSBarry Smith 1649a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 16509b94acceSBarry Smith 165131615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 165231615d04SLois Curfman McInnes SNESGetGradient() 16539b94acceSBarry Smith @*/ 16549b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 16559b94acceSBarry Smith { 165677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16579b94acceSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1, 165848d91487SBarry Smith "SNESGetFunction:For SNES_NONLINEAR_EQUATIONS only"); 16599b94acceSBarry Smith *r = snes->vec_func_always; 16609b94acceSBarry Smith return 0; 16619b94acceSBarry Smith } 16629b94acceSBarry Smith 16639b94acceSBarry Smith /*@C 16643638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 16659b94acceSBarry Smith 16669b94acceSBarry Smith Input Parameter: 16679b94acceSBarry Smith . snes - the SNES context 16689b94acceSBarry Smith 16699b94acceSBarry Smith Output Parameter: 16709b94acceSBarry Smith . r - the gradient 16719b94acceSBarry Smith 16729b94acceSBarry Smith Notes: 16739b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 16749b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 16759b94acceSBarry Smith SNESGetFunction(). 16769b94acceSBarry Smith 16779b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 16789b94acceSBarry Smith 167931615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 16809b94acceSBarry Smith @*/ 16819b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 16829b94acceSBarry Smith { 168377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16849b94acceSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 168548d91487SBarry Smith "SNESGetGradient:For SNES_UNCONSTRAINED_MINIMIZATION only"); 16869b94acceSBarry Smith *r = snes->vec_func_always; 16879b94acceSBarry Smith return 0; 16889b94acceSBarry Smith } 16899b94acceSBarry Smith 16909b94acceSBarry Smith /*@ 16919b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 16929b94acceSBarry Smith unconstrained minimization problems. 16939b94acceSBarry Smith 16949b94acceSBarry Smith Input Parameter: 16959b94acceSBarry Smith . snes - the SNES context 16969b94acceSBarry Smith 16979b94acceSBarry Smith Output Parameter: 16989b94acceSBarry Smith . r - the function 16999b94acceSBarry Smith 17009b94acceSBarry Smith Notes: 17019b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 17029b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 17039b94acceSBarry Smith SNESGetFunction(). 17049b94acceSBarry Smith 17059b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 17069b94acceSBarry Smith 170731615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 17089b94acceSBarry Smith @*/ 17099b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 17109b94acceSBarry Smith { 171177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 171274679c65SBarry Smith PetscValidScalarPointer(r); 17139b94acceSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1, 171448d91487SBarry Smith "SNESGetMinimizationFunction:For SNES_UNCONSTRAINED_MINIMIZATION only"); 17159b94acceSBarry Smith *r = snes->fc; 17169b94acceSBarry Smith return 0; 17179b94acceSBarry Smith } 17189b94acceSBarry Smith 17193c7409f5SSatish Balay /*@C 17203c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1721639f9d9dSBarry Smith SNES options in the database. You must NOT include the - at the beginning of 1722639f9d9dSBarry Smith the prefix name. 17233c7409f5SSatish Balay 17243c7409f5SSatish Balay Input Parameter: 17253c7409f5SSatish Balay . snes - the SNES context 17263c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 17273c7409f5SSatish Balay 17283c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1729a86d99e1SLois Curfman McInnes 1730a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 17313c7409f5SSatish Balay @*/ 17323c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 17333c7409f5SSatish Balay { 17343c7409f5SSatish Balay int ierr; 17353c7409f5SSatish Balay 173677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1737639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 17383c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 17393c7409f5SSatish Balay return 0; 17403c7409f5SSatish Balay } 17413c7409f5SSatish Balay 17423c7409f5SSatish Balay /*@C 1743f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1744639f9d9dSBarry Smith SNES options in the database. You must NOT include the - at the beginning of 1745639f9d9dSBarry Smith the prefix name. 17463c7409f5SSatish Balay 17473c7409f5SSatish Balay Input Parameter: 17483c7409f5SSatish Balay . snes - the SNES context 17493c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 17503c7409f5SSatish Balay 17513c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1752a86d99e1SLois Curfman McInnes 1753a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 17543c7409f5SSatish Balay @*/ 17553c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 17563c7409f5SSatish Balay { 17573c7409f5SSatish Balay int ierr; 17583c7409f5SSatish Balay 175977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1760639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 17613c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 17623c7409f5SSatish Balay return 0; 17633c7409f5SSatish Balay } 17643c7409f5SSatish Balay 1765c83590e2SSatish Balay /*@ 17663c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 17673c7409f5SSatish Balay SNES options in the database. 17683c7409f5SSatish Balay 17693c7409f5SSatish Balay Input Parameter: 17703c7409f5SSatish Balay . snes - the SNES context 17713c7409f5SSatish Balay 17723c7409f5SSatish Balay Output Parameter: 17733c7409f5SSatish Balay . prefix - pointer to the prefix string used 17743c7409f5SSatish Balay 17753c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 1776a86d99e1SLois Curfman McInnes 1777a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 17783c7409f5SSatish Balay @*/ 17793c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 17803c7409f5SSatish Balay { 17813c7409f5SSatish Balay int ierr; 17823c7409f5SSatish Balay 178377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1784639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 17853c7409f5SSatish Balay return 0; 17863c7409f5SSatish Balay } 17873c7409f5SSatish Balay 17883c7409f5SSatish Balay 17899b94acceSBarry Smith 17909b94acceSBarry Smith 17919b94acceSBarry Smith 1792