19b94acceSBarry Smith #ifndef lint 2*5eea60f9SBarry Smith static char vcid[] = "$Id: snes.c,v 1.116 1997/02/04 21:26:13 bsmith 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 1484cb2905SBarry Smith int SNESRegisterAllCalled = 0; 1584cb2905SBarry Smith 165615d1e5SSatish Balay #undef __FUNC__ 17*5eea60f9SBarry Smith #define __FUNC__ "SNESView" /* ADIC Ignore */ 189b94acceSBarry Smith /*@ 199b94acceSBarry Smith SNESView - Prints the SNES data structure. 209b94acceSBarry Smith 219b94acceSBarry Smith Input Parameters: 229b94acceSBarry Smith . SNES - the SNES context 239b94acceSBarry Smith . viewer - visualization context 249b94acceSBarry Smith 259b94acceSBarry Smith Options Database Key: 269b94acceSBarry Smith $ -snes_view : calls SNESView() at end of SNESSolve() 279b94acceSBarry Smith 289b94acceSBarry Smith Notes: 299b94acceSBarry Smith The available visualization contexts include 306d4a8577SBarry Smith $ VIEWER_STDOUT_SELF - standard output (default) 316d4a8577SBarry Smith $ VIEWER_STDOUT_WORLD - synchronized standard 329b94acceSBarry Smith $ output where only the first processor opens 339b94acceSBarry Smith $ the file. All other processors send their 349b94acceSBarry Smith $ data to the first processor to print. 359b94acceSBarry Smith 369b94acceSBarry Smith The user can open alternative vistualization contexts with 37dbb450caSBarry Smith $ ViewerFileOpenASCII() - output to a specified file 389b94acceSBarry Smith 399b94acceSBarry Smith .keywords: SNES, view 409b94acceSBarry Smith 41dbb450caSBarry Smith .seealso: ViewerFileOpenASCII() 429b94acceSBarry Smith @*/ 439b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 449b94acceSBarry Smith { 459b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 469b94acceSBarry Smith FILE *fd; 479b94acceSBarry Smith int ierr; 489b94acceSBarry Smith SLES sles; 499b94acceSBarry Smith char *method; 5019bcc07fSBarry Smith ViewerType vtype; 519b94acceSBarry Smith 5274679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5374679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5474679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5574679c65SBarry Smith 5619bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 5719bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 5890ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 5977c4ece6SBarry Smith PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 604b0e389bSBarry Smith SNESGetType(snes,PETSC_NULL,&method); 6177c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," method: %s\n",method); 629b94acceSBarry Smith if (snes->view) (*snes->view)((PetscObject)snes,viewer); 6377c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 649b94acceSBarry Smith " maximum iterations=%d, maximum function evaluations=%d\n", 659b94acceSBarry Smith snes->max_its,snes->max_funcs); 6677c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 679b94acceSBarry Smith " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 689b94acceSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol); 697a00f4a9SLois Curfman McInnes PetscFPrintf(snes->comm,fd, 707a00f4a9SLois Curfman McInnes " total number of linear solver iterations=%d\n",snes->linear_its); 71ae3c334cSLois Curfman McInnes PetscFPrintf(snes->comm,fd, 7268632166SLois Curfman McInnes " total number of function evaluations=%d\n",snes->nfuncs); 739b94acceSBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 7477c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 759b94acceSBarry Smith if (snes->ksp_ewconv) { 769b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 779b94acceSBarry Smith if (kctx) { 7877c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 799b94acceSBarry Smith " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 809b94acceSBarry Smith kctx->version); 8177c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 829b94acceSBarry Smith " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 839b94acceSBarry Smith kctx->rtol_max,kctx->threshold); 8477c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 859b94acceSBarry Smith kctx->gamma,kctx->alpha,kctx->alpha2); 869b94acceSBarry Smith } 879b94acceSBarry Smith } 88c30f285eSLois Curfman McInnes } else if (vtype == STRING_VIEWER) { 89c30f285eSLois Curfman McInnes SNESGetType(snes,PETSC_NULL,&method); 90c30f285eSLois Curfman McInnes ViewerStringSPrintf(viewer," %-3.3s",method); 9119bcc07fSBarry Smith } 929b94acceSBarry Smith SNESGetSLES(snes,&sles); 939b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 949b94acceSBarry Smith return 0; 959b94acceSBarry Smith } 969b94acceSBarry Smith 97639f9d9dSBarry Smith /* 98639f9d9dSBarry Smith We retain a list of functions that also take SNES command 99639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 100639f9d9dSBarry Smith */ 101639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 102639f9d9dSBarry Smith static int numberofsetfromoptions; 103639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 104639f9d9dSBarry Smith 1055615d1e5SSatish Balay #undef __FUNC__ 106*5eea60f9SBarry Smith #define __FUNC__ "SNESAddOptionsChecker" /* ADIC Ignore */ 107639f9d9dSBarry Smith /*@ 108639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 109639f9d9dSBarry Smith 110639f9d9dSBarry Smith Input Parameter: 111639f9d9dSBarry Smith . snescheck - function that checks for options 112639f9d9dSBarry Smith 113639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 114639f9d9dSBarry Smith @*/ 115639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 116639f9d9dSBarry Smith { 117639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 118e3372554SBarry Smith SETERRQ(1,0,"Too many options checkers, only 5 allowed"); 119639f9d9dSBarry Smith } 120639f9d9dSBarry Smith 121639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 122639f9d9dSBarry Smith return 0; 123639f9d9dSBarry Smith } 124639f9d9dSBarry Smith 1255615d1e5SSatish Balay #undef __FUNC__ 1265615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1279b94acceSBarry Smith /*@ 128682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1299b94acceSBarry Smith 1309b94acceSBarry Smith Input Parameter: 1319b94acceSBarry Smith . snes - the SNES context 1329b94acceSBarry Smith 1339b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1349b94acceSBarry Smith 135a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1369b94acceSBarry Smith @*/ 1379b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1389b94acceSBarry Smith { 1394b0e389bSBarry Smith SNESType method; 1409b94acceSBarry Smith double tmp; 1419b94acceSBarry Smith SLES sles; 14281f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1433c7409f5SSatish Balay int version = PETSC_DEFAULT; 1449b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1459b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1469b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1479b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1489b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1499b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1509b94acceSBarry Smith 15181f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 15281f57ec7SBarry Smith 15381f57ec7SBarry Smith 15477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 155e3372554SBarry Smith if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp"); 156052efed2SBarry Smith ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 157052efed2SBarry Smith if (flg) { 158052efed2SBarry Smith ierr = SNESSetType(snes,method); CHKERRQ(ierr); 1599b94acceSBarry Smith } 1605334005bSBarry Smith else if (!snes->set_method_called) { 1615334005bSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 16240191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 1635334005bSBarry Smith } 1645334005bSBarry Smith else { 16540191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 1665334005bSBarry Smith } 1675334005bSBarry Smith } 1685334005bSBarry Smith 1693c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1703c7409f5SSatish Balay if (flg) { SNESPrintHelp(snes); } 1713c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 172d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);} 1733c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 174d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);} 1753c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 176d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); } 1773c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 1783c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 179d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 180d7a720efSLois Curfman McInnes if (flg) { SNESSetTrustRegionTolerance(snes,tmp); } 181d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 182d7a720efSLois Curfman McInnes if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); } 1833c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 1843c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 185b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 186b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 187b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 188b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 189b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 190b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 191b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 1929b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 1939b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 1945bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 1955bbad29bSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);} 1963c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 197639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 1983c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 199d31a3109SLois Curfman McInnes if (flg) { 200639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 201d31a3109SLois Curfman McInnes } 20281f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2033c7409f5SSatish Balay if (flg) { 20417699dbbSLois Curfman McInnes int rank = 0; 205d7e8b826SBarry Smith DrawLG lg; 20617699dbbSLois Curfman McInnes MPI_Initialized(&rank); 20717699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 20817699dbbSLois Curfman McInnes if (!rank) { 20981f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2109b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 211f8c826e1SBarry Smith snes->xmonitor = lg; 2129b94acceSBarry Smith PLogObjectParent(snes,lg); 2139b94acceSBarry Smith } 2149b94acceSBarry Smith } 2153c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2163c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2179b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2189b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 21994a424c1SBarry Smith PLogInfo(snes, 22031615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 22131615d04SLois Curfman McInnes } 22231615d04SLois Curfman McInnes else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 22331615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 22431615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 22594a424c1SBarry Smith PLogInfo(snes, 22631615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2279b94acceSBarry Smith } 228639f9d9dSBarry Smith 229639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 230639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 231639f9d9dSBarry Smith } 232639f9d9dSBarry Smith 2339b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2349b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 2359b94acceSBarry Smith if (!snes->setfromoptions) return 0; 2369b94acceSBarry Smith return (*snes->setfromoptions)(snes); 2379b94acceSBarry Smith } 2389b94acceSBarry Smith 2395615d1e5SSatish Balay #undef __FUNC__ 240*5eea60f9SBarry Smith #define __FUNC__ "SNESPrintHelp" /* ADIC Ignore */ 2419b94acceSBarry Smith /*@ 2429b94acceSBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2439b94acceSBarry Smith 2449b94acceSBarry Smith Input Parameter: 2459b94acceSBarry Smith . snes - the SNES context 2469b94acceSBarry Smith 247a703fe33SLois Curfman McInnes Options Database Keys: 248a703fe33SLois Curfman McInnes $ -help, -h 249a703fe33SLois Curfman McInnes 2509b94acceSBarry Smith .keywords: SNES, nonlinear, help 2519b94acceSBarry Smith 252682d7d0cSBarry Smith .seealso: SNESSetFromOptions() 2539b94acceSBarry Smith @*/ 2549b94acceSBarry Smith int SNESPrintHelp(SNES snes) 2559b94acceSBarry Smith { 2566daaf66cSBarry Smith char p[64]; 2576daaf66cSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2586daaf66cSBarry Smith 25977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2606daaf66cSBarry Smith 2616daaf66cSBarry Smith PetscStrcpy(p,"-"); 2626daaf66cSBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2636daaf66cSBarry Smith 2646daaf66cSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2656daaf66cSBarry Smith 266d31a3109SLois Curfman McInnes PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 26733455573SSatish Balay SNESPrintTypes_Private(snes->comm,p,"snes_type"); 26877c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 269b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 270b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 271b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 272b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 273b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 274b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 2755bbad29bSBarry Smith PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 2765bbad29bSBarry Smith PetscPrintf(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 277d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor\n",p); 278d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 279b18e04deSLois Curfman McInnes if (snes->type == SNES_NONLINEAR_EQUATIONS) { 28077c4ece6SBarry Smith PetscPrintf(snes->comm, 281d31a3109SLois Curfman McInnes " Options for solving systems of nonlinear equations only:\n"); 28277c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 28377c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 284ef1dfb11SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 28577c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 28677c4ece6SBarry Smith PetscPrintf(snes->comm, 287b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 28877c4ece6SBarry Smith PetscPrintf(snes->comm, 289b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 29077c4ece6SBarry Smith PetscPrintf(snes->comm, 291b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 29277c4ece6SBarry Smith PetscPrintf(snes->comm, 293b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 29477c4ece6SBarry Smith PetscPrintf(snes->comm, 295b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 29677c4ece6SBarry Smith PetscPrintf(snes->comm, 297b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 29877c4ece6SBarry Smith PetscPrintf(snes->comm, 299b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 300b18e04deSLois Curfman McInnes } 301b18e04deSLois Curfman McInnes else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 30277c4ece6SBarry Smith PetscPrintf(snes->comm, 303d31a3109SLois Curfman McInnes " Options for solving unconstrained minimization problems only:\n"); 304b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 30577c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 30677c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 307b18e04deSLois Curfman McInnes } 3084537a946SLois Curfman McInnes PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 30977c4ece6SBarry Smith PetscPrintf(snes->comm,"a particular method\n"); 3106daaf66cSBarry Smith if (snes->printhelp) (*snes->printhelp)(snes,p); 3119b94acceSBarry Smith return 0; 3129b94acceSBarry Smith } 313a847f771SSatish Balay 3145615d1e5SSatish Balay #undef __FUNC__ 315*5eea60f9SBarry Smith #define __FUNC__ "SNESSetApplicationContext" /* ADIC Ignore */ 3169b94acceSBarry Smith /*@ 3179b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 3189b94acceSBarry Smith the nonlinear solvers. 3199b94acceSBarry Smith 3209b94acceSBarry Smith Input Parameters: 3219b94acceSBarry Smith . snes - the SNES context 3229b94acceSBarry Smith . usrP - optional user context 3239b94acceSBarry Smith 3249b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3259b94acceSBarry Smith 3269b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3279b94acceSBarry Smith @*/ 3289b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3299b94acceSBarry Smith { 33077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3319b94acceSBarry Smith snes->user = usrP; 3329b94acceSBarry Smith return 0; 3339b94acceSBarry Smith } 33474679c65SBarry Smith 3355615d1e5SSatish Balay #undef __FUNC__ 336*5eea60f9SBarry Smith #define __FUNC__ "SNESGetApplicationContext" /* ADIC Ignore */ 3379b94acceSBarry Smith /*@C 3389b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3399b94acceSBarry Smith nonlinear solvers. 3409b94acceSBarry Smith 3419b94acceSBarry Smith Input Parameter: 3429b94acceSBarry Smith . snes - SNES context 3439b94acceSBarry Smith 3449b94acceSBarry Smith Output Parameter: 3459b94acceSBarry Smith . usrP - user context 3469b94acceSBarry Smith 3479b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3489b94acceSBarry Smith 3499b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3509b94acceSBarry Smith @*/ 3519b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3529b94acceSBarry Smith { 35377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3549b94acceSBarry Smith *usrP = snes->user; 3559b94acceSBarry Smith return 0; 3569b94acceSBarry Smith } 35774679c65SBarry Smith 3585615d1e5SSatish Balay #undef __FUNC__ 359*5eea60f9SBarry Smith #define __FUNC__ "SNESGetIterationNumber" /* ADIC Ignore */ 3609b94acceSBarry Smith /*@ 3619b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3629b94acceSBarry Smith nonlinear solver. 3639b94acceSBarry Smith 3649b94acceSBarry Smith Input Parameter: 3659b94acceSBarry Smith . snes - SNES context 3669b94acceSBarry Smith 3679b94acceSBarry Smith Output Parameter: 3689b94acceSBarry Smith . iter - iteration number 3699b94acceSBarry Smith 3709b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3719b94acceSBarry Smith @*/ 3729b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3739b94acceSBarry Smith { 37477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 37574679c65SBarry Smith PetscValidIntPointer(iter); 3769b94acceSBarry Smith *iter = snes->iter; 3779b94acceSBarry Smith return 0; 3789b94acceSBarry Smith } 37974679c65SBarry Smith 3805615d1e5SSatish Balay #undef __FUNC__ 3815615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 3829b94acceSBarry Smith /*@ 3839b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3849b94acceSBarry Smith with SNESSSetFunction(). 3859b94acceSBarry Smith 3869b94acceSBarry Smith Input Parameter: 3879b94acceSBarry Smith . snes - SNES context 3889b94acceSBarry Smith 3899b94acceSBarry Smith Output Parameter: 3909b94acceSBarry Smith . fnorm - 2-norm of function 3919b94acceSBarry Smith 3929b94acceSBarry Smith Note: 3939b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 394a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 395a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 3969b94acceSBarry Smith 3979b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 398a86d99e1SLois Curfman McInnes 399a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 4009b94acceSBarry Smith @*/ 4019b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 4029b94acceSBarry Smith { 40377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 40474679c65SBarry Smith PetscValidScalarPointer(fnorm); 40574679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 406d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 40774679c65SBarry Smith } 4089b94acceSBarry Smith *fnorm = snes->norm; 4099b94acceSBarry Smith return 0; 4109b94acceSBarry Smith } 41174679c65SBarry Smith 4125615d1e5SSatish Balay #undef __FUNC__ 4135615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 4149b94acceSBarry Smith /*@ 4159b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 4169b94acceSBarry Smith with SNESSSetGradient(). 4179b94acceSBarry Smith 4189b94acceSBarry Smith Input Parameter: 4199b94acceSBarry Smith . snes - SNES context 4209b94acceSBarry Smith 4219b94acceSBarry Smith Output Parameter: 4229b94acceSBarry Smith . fnorm - 2-norm of gradient 4239b94acceSBarry Smith 4249b94acceSBarry Smith Note: 4259b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 426a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 427a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4289b94acceSBarry Smith 4299b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 430a86d99e1SLois Curfman McInnes 431a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4329b94acceSBarry Smith @*/ 4339b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4349b94acceSBarry Smith { 43577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 43674679c65SBarry Smith PetscValidScalarPointer(gnorm); 43774679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 438d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 43974679c65SBarry Smith } 4409b94acceSBarry Smith *gnorm = snes->norm; 4419b94acceSBarry Smith return 0; 4429b94acceSBarry Smith } 44374679c65SBarry Smith 4445615d1e5SSatish Balay #undef __FUNC__ 445*5eea60f9SBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" /* ADIC Ignore */ 4469b94acceSBarry Smith /*@ 4479b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4489b94acceSBarry Smith attempted by the nonlinear solver. 4499b94acceSBarry Smith 4509b94acceSBarry Smith Input Parameter: 4519b94acceSBarry Smith . snes - SNES context 4529b94acceSBarry Smith 4539b94acceSBarry Smith Output Parameter: 4549b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4559b94acceSBarry Smith 456c96a6f78SLois Curfman McInnes Notes: 457c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 458c96a6f78SLois Curfman McInnes 4599b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4609b94acceSBarry Smith @*/ 4619b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4629b94acceSBarry Smith { 46377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 46474679c65SBarry Smith PetscValidIntPointer(nfails); 4659b94acceSBarry Smith *nfails = snes->nfailures; 4669b94acceSBarry Smith return 0; 4679b94acceSBarry Smith } 468a847f771SSatish Balay 4695615d1e5SSatish Balay #undef __FUNC__ 470*5eea60f9SBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" /* ADIC Ignore */ 471c96a6f78SLois Curfman McInnes /*@ 472c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 473c96a6f78SLois Curfman McInnes used by the nonlinear solver. 474c96a6f78SLois Curfman McInnes 475c96a6f78SLois Curfman McInnes Input Parameter: 476c96a6f78SLois Curfman McInnes . snes - SNES context 477c96a6f78SLois Curfman McInnes 478c96a6f78SLois Curfman McInnes Output Parameter: 479c96a6f78SLois Curfman McInnes . lits - number of linear iterations 480c96a6f78SLois Curfman McInnes 481c96a6f78SLois Curfman McInnes Notes: 482c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 483c96a6f78SLois Curfman McInnes 484c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 485c96a6f78SLois Curfman McInnes @*/ 48686bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 487c96a6f78SLois Curfman McInnes { 488c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 489c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 490c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 491c96a6f78SLois Curfman McInnes return 0; 492c96a6f78SLois Curfman McInnes } 493c96a6f78SLois Curfman McInnes 494c96a6f78SLois Curfman McInnes #undef __FUNC__ 495*5eea60f9SBarry Smith #define __FUNC__ "SNESGetSLES" /* ADIC Ignore */ 4969b94acceSBarry Smith /*@C 4979b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 4989b94acceSBarry Smith 4999b94acceSBarry Smith Input Parameter: 5009b94acceSBarry Smith . snes - the SNES context 5019b94acceSBarry Smith 5029b94acceSBarry Smith Output Parameter: 5039b94acceSBarry Smith . sles - the SLES context 5049b94acceSBarry Smith 5059b94acceSBarry Smith Notes: 5069b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 5079b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5089b94acceSBarry Smith KSP and PC contexts as well. 5099b94acceSBarry Smith 5109b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 5119b94acceSBarry Smith 5129b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 5139b94acceSBarry Smith @*/ 5149b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 5159b94acceSBarry Smith { 51677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5179b94acceSBarry Smith *sles = snes->sles; 5189b94acceSBarry Smith return 0; 5199b94acceSBarry Smith } 5209b94acceSBarry Smith /* -----------------------------------------------------------*/ 5215615d1e5SSatish Balay #undef __FUNC__ 5225615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 5239b94acceSBarry Smith /*@C 5249b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5259b94acceSBarry Smith 5269b94acceSBarry Smith Input Parameter: 5279b94acceSBarry Smith . comm - MPI communicator 5289b94acceSBarry Smith . type - type of method, one of 5299b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS 5309b94acceSBarry Smith $ (for systems of nonlinear equations) 5319b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION 5329b94acceSBarry Smith $ (for unconstrained minimization) 5339b94acceSBarry Smith 5349b94acceSBarry Smith Output Parameter: 5359b94acceSBarry Smith . outsnes - the new SNES context 5369b94acceSBarry Smith 537c1f60f51SBarry Smith Options Database Key: 53819bd6747SLois Curfman McInnes $ -snes_mf - use default matrix-free Jacobian-vector products, 53919bd6747SLois Curfman McInnes $ and no preconditioning matrix 54019bd6747SLois Curfman McInnes $ -snes_mf_operator - use default matrix-free Jacobian-vector 54119bd6747SLois Curfman McInnes $ products, and a user-provided preconditioning matrix 54219bd6747SLois Curfman McInnes $ as set by SNESSetJacobian() 54319bd6747SLois Curfman McInnes $ -snes_fd - use (slow!) finite differences to compute Jacobian 544c1f60f51SBarry Smith 5459b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5469b94acceSBarry Smith 54763a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 5489b94acceSBarry Smith @*/ 5494b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 5509b94acceSBarry Smith { 5519b94acceSBarry Smith int ierr; 5529b94acceSBarry Smith SNES snes; 5539b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 55437fcc0dbSBarry Smith 5559b94acceSBarry Smith *outsnes = 0; 5562a0acf01SLois Curfman McInnes if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 557d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 5580452661fSBarry Smith PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm); 5599b94acceSBarry Smith PLogObjectCreate(snes); 5609b94acceSBarry Smith snes->max_its = 50; 5619b94acceSBarry Smith snes->max_funcs = 1000; 5629b94acceSBarry Smith snes->norm = 0.0; 5635a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 564b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 565b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5669b94acceSBarry Smith snes->atol = 1.e-10; 5675a2d0531SBarry Smith } 568b4874afaSBarry Smith else { 569b4874afaSBarry Smith snes->rtol = 1.e-8; 570b4874afaSBarry Smith snes->ttol = 0.0; 571b4874afaSBarry Smith snes->atol = 1.e-50; 572b4874afaSBarry Smith } 5739b94acceSBarry Smith snes->xtol = 1.e-8; 574b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5759b94acceSBarry Smith snes->nfuncs = 0; 5769b94acceSBarry Smith snes->nfailures = 0; 5777a00f4a9SLois Curfman McInnes snes->linear_its = 0; 578639f9d9dSBarry Smith snes->numbermonitors = 0; 5799b94acceSBarry Smith snes->data = 0; 5809b94acceSBarry Smith snes->view = 0; 5819b94acceSBarry Smith snes->computeumfunction = 0; 5829b94acceSBarry Smith snes->umfunP = 0; 5839b94acceSBarry Smith snes->fc = 0; 5849b94acceSBarry Smith snes->deltatol = 1.e-12; 5859b94acceSBarry Smith snes->fmin = -1.e30; 5869b94acceSBarry Smith snes->method_class = type; 5879b94acceSBarry Smith snes->set_method_called = 0; 588a703fe33SLois Curfman McInnes snes->setup_called = 0; 5899b94acceSBarry Smith snes->ksp_ewconv = 0; 5907d1a2b2bSBarry Smith snes->type = -1; 5916f24a144SLois Curfman McInnes snes->xmonitor = 0; 5926f24a144SLois Curfman McInnes snes->vwork = 0; 5936f24a144SLois Curfman McInnes snes->nwork = 0; 594c9005455SLois Curfman McInnes snes->conv_hist_size = 0; 595c9005455SLois Curfman McInnes snes->conv_act_size = 0; 596c9005455SLois Curfman McInnes snes->conv_hist = 0; 5979b94acceSBarry Smith 5989b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 5990452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 6009b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6019b94acceSBarry Smith kctx->version = 2; 6029b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6039b94acceSBarry Smith this was too large for some test cases */ 6049b94acceSBarry Smith kctx->rtol_last = 0; 6059b94acceSBarry Smith kctx->rtol_max = .9; 6069b94acceSBarry Smith kctx->gamma = 1.0; 6079b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6089b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6099b94acceSBarry Smith kctx->threshold = .1; 6109b94acceSBarry Smith kctx->lresid_last = 0; 6119b94acceSBarry Smith kctx->norm_last = 0; 6129b94acceSBarry Smith 6139b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 6149b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 6155334005bSBarry Smith 6169b94acceSBarry Smith *outsnes = snes; 6179b94acceSBarry Smith return 0; 6189b94acceSBarry Smith } 6199b94acceSBarry Smith 6209b94acceSBarry Smith /* --------------------------------------------------------------- */ 6215615d1e5SSatish Balay #undef __FUNC__ 622*5eea60f9SBarry Smith #define __FUNC__ "SNESSetFunction" /* ADIC Ignore */ 6239b94acceSBarry Smith /*@C 6249b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6259b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6269b94acceSBarry Smith equations. 6279b94acceSBarry Smith 6289b94acceSBarry Smith Input Parameters: 6299b94acceSBarry Smith . snes - the SNES context 6309b94acceSBarry Smith . func - function evaluation routine 6319b94acceSBarry Smith . r - vector to store function value 6322cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 6332cd2dfdcSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6349b94acceSBarry Smith 6359b94acceSBarry Smith Calling sequence of func: 636313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Vec f,void *ctx); 6379b94acceSBarry Smith 6389b94acceSBarry Smith . x - input vector 639313e4042SLois Curfman McInnes . f - function vector 6402cd2dfdcSLois Curfman McInnes . ctx - optional user-defined function context 6419b94acceSBarry Smith 6429b94acceSBarry Smith Notes: 6439b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6449b94acceSBarry Smith $ f'(x) x = -f(x), 6459b94acceSBarry Smith $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 6469b94acceSBarry Smith 6479b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6489b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6499b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 6509b94acceSBarry Smith 6519b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6529b94acceSBarry Smith 653a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6549b94acceSBarry Smith @*/ 6555334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 6569b94acceSBarry Smith { 65777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 658e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6599b94acceSBarry Smith snes->computefunction = func; 6609b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6619b94acceSBarry Smith snes->funP = ctx; 6629b94acceSBarry Smith return 0; 6639b94acceSBarry Smith } 6649b94acceSBarry Smith 6655615d1e5SSatish Balay #undef __FUNC__ 6665615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 6679b94acceSBarry Smith /*@ 6689b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6699b94acceSBarry Smith SNESSetFunction(). 6709b94acceSBarry Smith 6719b94acceSBarry Smith Input Parameters: 6729b94acceSBarry Smith . snes - the SNES context 6739b94acceSBarry Smith . x - input vector 6749b94acceSBarry Smith 6759b94acceSBarry Smith Output Parameter: 6763638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6779b94acceSBarry Smith 6781bffabb2SLois Curfman McInnes Notes: 6799b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6809b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6819b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6829b94acceSBarry Smith 6839b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6849b94acceSBarry Smith 685a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6869b94acceSBarry Smith @*/ 6879b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 6889b94acceSBarry Smith { 6899b94acceSBarry Smith int ierr; 6909b94acceSBarry Smith 69174679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 692e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6939b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 6949b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 695ae3c334cSLois Curfman McInnes snes->nfuncs++; 6969b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 6979b94acceSBarry Smith return 0; 6989b94acceSBarry Smith } 6999b94acceSBarry Smith 7005615d1e5SSatish Balay #undef __FUNC__ 701*5eea60f9SBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" /* ADIC Ignore */ 7029b94acceSBarry Smith /*@C 7039b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 7049b94acceSBarry Smith unconstrained minimization. 7059b94acceSBarry Smith 7069b94acceSBarry Smith Input Parameters: 7079b94acceSBarry Smith . snes - the SNES context 7089b94acceSBarry Smith . func - function evaluation routine 709044dda88SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 710044dda88SLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 7119b94acceSBarry Smith 7129b94acceSBarry Smith Calling sequence of func: 7139b94acceSBarry Smith . func (SNES snes,Vec x,double *f,void *ctx); 7149b94acceSBarry Smith 7159b94acceSBarry Smith . x - input vector 7169b94acceSBarry Smith . f - function 717044dda88SLois Curfman McInnes . ctx - [optional] user-defined function context 7189b94acceSBarry Smith 7199b94acceSBarry Smith Notes: 7209b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7219b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7229b94acceSBarry Smith SNESSetFunction(). 7239b94acceSBarry Smith 7249b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 7259b94acceSBarry Smith 726a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 727a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 7289b94acceSBarry Smith @*/ 7299b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 7309b94acceSBarry Smith void *ctx) 7319b94acceSBarry Smith { 73277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 733e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7349b94acceSBarry Smith snes->computeumfunction = func; 7359b94acceSBarry Smith snes->umfunP = ctx; 7369b94acceSBarry Smith return 0; 7379b94acceSBarry Smith } 7389b94acceSBarry Smith 7395615d1e5SSatish Balay #undef __FUNC__ 7405615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 7419b94acceSBarry Smith /*@ 7429b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 7439b94acceSBarry Smith set with SNESSetMinimizationFunction(). 7449b94acceSBarry Smith 7459b94acceSBarry Smith Input Parameters: 7469b94acceSBarry Smith . snes - the SNES context 7479b94acceSBarry Smith . x - input vector 7489b94acceSBarry Smith 7499b94acceSBarry Smith Output Parameter: 7509b94acceSBarry Smith . y - function value 7519b94acceSBarry Smith 7529b94acceSBarry Smith Notes: 7539b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 7549b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7559b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 756a86d99e1SLois Curfman McInnes 757a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 758a86d99e1SLois Curfman McInnes 759a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 760a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 7619b94acceSBarry Smith @*/ 7629b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 7639b94acceSBarry Smith { 7649b94acceSBarry Smith int ierr; 765e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7669b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 7679b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 768ae3c334cSLois Curfman McInnes snes->nfuncs++; 7699b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7709b94acceSBarry Smith return 0; 7719b94acceSBarry Smith } 7729b94acceSBarry Smith 7735615d1e5SSatish Balay #undef __FUNC__ 774*5eea60f9SBarry Smith #define __FUNC__ "SNESSetGradient" /* ADIC Ignore */ 7759b94acceSBarry Smith /*@C 7769b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 7779b94acceSBarry Smith vector for use by the SNES routines. 7789b94acceSBarry Smith 7799b94acceSBarry Smith Input Parameters: 7809b94acceSBarry Smith . snes - the SNES context 7819b94acceSBarry Smith . func - function evaluation routine 782044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 783044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 7849b94acceSBarry Smith . r - vector to store gradient value 7859b94acceSBarry Smith 7869b94acceSBarry Smith Calling sequence of func: 7879b94acceSBarry Smith . func (SNES, Vec x, Vec g, void *ctx); 7889b94acceSBarry Smith 7899b94acceSBarry Smith . x - input vector 7909b94acceSBarry Smith . g - gradient vector 791044dda88SLois Curfman McInnes . ctx - optional user-defined gradient context 7929b94acceSBarry Smith 7939b94acceSBarry Smith Notes: 7949b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7959b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7969b94acceSBarry Smith SNESSetFunction(). 7979b94acceSBarry Smith 7989b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 7999b94acceSBarry Smith 800a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 801a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 8029b94acceSBarry Smith @*/ 80374679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 8049b94acceSBarry Smith { 80577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 806e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8079b94acceSBarry Smith snes->computefunction = func; 8089b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 8099b94acceSBarry Smith snes->funP = ctx; 8109b94acceSBarry Smith return 0; 8119b94acceSBarry Smith } 8129b94acceSBarry Smith 8135615d1e5SSatish Balay #undef __FUNC__ 8145615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 8159b94acceSBarry Smith /*@ 816a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 817a86d99e1SLois Curfman McInnes SNESSetGradient(). 8189b94acceSBarry Smith 8199b94acceSBarry Smith Input Parameters: 8209b94acceSBarry Smith . snes - the SNES context 8219b94acceSBarry Smith . x - input vector 8229b94acceSBarry Smith 8239b94acceSBarry Smith Output Parameter: 8249b94acceSBarry Smith . y - gradient vector 8259b94acceSBarry Smith 8269b94acceSBarry Smith Notes: 8279b94acceSBarry Smith SNESComputeGradient() is valid only for 8289b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 8299b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 830a86d99e1SLois Curfman McInnes 831a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 832a86d99e1SLois Curfman McInnes 833a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 834a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 8359b94acceSBarry Smith @*/ 8369b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 8379b94acceSBarry Smith { 8389b94acceSBarry Smith int ierr; 83974679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 840e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8419b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 8429b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 8439b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 8449b94acceSBarry Smith return 0; 8459b94acceSBarry Smith } 8469b94acceSBarry Smith 8475615d1e5SSatish Balay #undef __FUNC__ 8485615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 84962fef451SLois Curfman McInnes /*@ 85062fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 85162fef451SLois Curfman McInnes set with SNESSetJacobian(). 85262fef451SLois Curfman McInnes 85362fef451SLois Curfman McInnes Input Parameters: 85462fef451SLois Curfman McInnes . snes - the SNES context 85562fef451SLois Curfman McInnes . x - input vector 85662fef451SLois Curfman McInnes 85762fef451SLois Curfman McInnes Output Parameters: 85862fef451SLois Curfman McInnes . A - Jacobian matrix 85962fef451SLois Curfman McInnes . B - optional preconditioning matrix 86062fef451SLois Curfman McInnes . flag - flag indicating matrix structure 86162fef451SLois Curfman McInnes 86262fef451SLois Curfman McInnes Notes: 86362fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 86462fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 86562fef451SLois Curfman McInnes 866dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 867dc5a77f8SLois Curfman McInnes flag parameter. 86862fef451SLois Curfman McInnes 86962fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 87062fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 87162fef451SLois Curfman McInnes methods is SNESComputeJacobian(). 87262fef451SLois Curfman McInnes 87362fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 87462fef451SLois Curfman McInnes 87562fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 87662fef451SLois Curfman McInnes @*/ 8779b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8789b94acceSBarry Smith { 8799b94acceSBarry Smith int ierr; 88074679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 881e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 8829b94acceSBarry Smith if (!snes->computejacobian) return 0; 8839b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 884c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 8859b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 8869b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 8876d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 88877c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 88977c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 8909b94acceSBarry Smith return 0; 8919b94acceSBarry Smith } 8929b94acceSBarry Smith 8935615d1e5SSatish Balay #undef __FUNC__ 8945615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 89562fef451SLois Curfman McInnes /*@ 89662fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 89762fef451SLois Curfman McInnes set with SNESSetHessian(). 89862fef451SLois Curfman McInnes 89962fef451SLois Curfman McInnes Input Parameters: 90062fef451SLois Curfman McInnes . snes - the SNES context 90162fef451SLois Curfman McInnes . x - input vector 90262fef451SLois Curfman McInnes 90362fef451SLois Curfman McInnes Output Parameters: 90462fef451SLois Curfman McInnes . A - Hessian matrix 90562fef451SLois Curfman McInnes . B - optional preconditioning matrix 90662fef451SLois Curfman McInnes . flag - flag indicating matrix structure 90762fef451SLois Curfman McInnes 90862fef451SLois Curfman McInnes Notes: 90962fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 91062fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 91162fef451SLois Curfman McInnes 912dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 913dc5a77f8SLois Curfman McInnes flag parameter. 91462fef451SLois Curfman McInnes 91562fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 91662fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 91762fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 91862fef451SLois Curfman McInnes 91962fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 92062fef451SLois Curfman McInnes 921a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 922a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 92362fef451SLois Curfman McInnes @*/ 92462fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 9259b94acceSBarry Smith { 9269b94acceSBarry Smith int ierr; 92774679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 928e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 9299b94acceSBarry Smith if (!snes->computejacobian) return 0; 93062fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 9310f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 93262fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 93362fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 9340f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 93577c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 93677c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9379b94acceSBarry Smith return 0; 9389b94acceSBarry Smith } 9399b94acceSBarry Smith 9405615d1e5SSatish Balay #undef __FUNC__ 941*5eea60f9SBarry Smith #define __FUNC__ "SNESSetJacobian" /* ADIC Ignore */ 9429b94acceSBarry Smith /*@C 9439b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 944044dda88SLois Curfman McInnes location to store the matrix. 9459b94acceSBarry Smith 9469b94acceSBarry Smith Input Parameters: 9479b94acceSBarry Smith . snes - the SNES context 9489b94acceSBarry Smith . A - Jacobian matrix 9499b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 9509b94acceSBarry Smith . func - Jacobian evaluation routine 9512cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 9522cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 9539b94acceSBarry Smith 9549b94acceSBarry Smith Calling sequence of func: 955313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9569b94acceSBarry Smith 9579b94acceSBarry Smith . x - input vector 9589b94acceSBarry Smith . A - Jacobian matrix 9599b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 960ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 961ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 9622cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Jacobian context 9639b94acceSBarry Smith 9649b94acceSBarry Smith Notes: 965dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 9662cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 967ac21db08SLois Curfman McInnes 968ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9699b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9709b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9719b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9729b94acceSBarry Smith throughout the global iterations. 9739b94acceSBarry Smith 9749b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9759b94acceSBarry Smith 976ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 9779b94acceSBarry Smith @*/ 9789b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 9799b94acceSBarry Smith MatStructure*,void*),void *ctx) 9809b94acceSBarry Smith { 98177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 98274679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 983e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 9849b94acceSBarry Smith snes->computejacobian = func; 9859b94acceSBarry Smith snes->jacP = ctx; 9869b94acceSBarry Smith snes->jacobian = A; 9879b94acceSBarry Smith snes->jacobian_pre = B; 9889b94acceSBarry Smith return 0; 9899b94acceSBarry Smith } 99062fef451SLois Curfman McInnes 9915615d1e5SSatish Balay #undef __FUNC__ 992*5eea60f9SBarry Smith #define __FUNC__ "SNESGetJacobian" /* ADIC Ignore */ 993b4fd4287SBarry Smith /*@ 994b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 995b4fd4287SBarry Smith provided context for evaluating the Jacobian. 996b4fd4287SBarry Smith 997b4fd4287SBarry Smith Input Parameter: 998b4fd4287SBarry Smith . snes - the nonlinear solver context 999b4fd4287SBarry Smith 1000b4fd4287SBarry Smith Output Parameters: 1001b4fd4287SBarry Smith . A - location to stash Jacobian matrix (or PETSC_NULL) 1002b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 1003b4fd4287SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 1004b4fd4287SBarry Smith 1005b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 1006b4fd4287SBarry Smith @*/ 1007b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1008b4fd4287SBarry Smith { 100974679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 1010e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 1011b4fd4287SBarry Smith if (A) *A = snes->jacobian; 1012b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 1013b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 1014b4fd4287SBarry Smith return 0; 1015b4fd4287SBarry Smith } 1016b4fd4287SBarry Smith 10175615d1e5SSatish Balay #undef __FUNC__ 1018*5eea60f9SBarry Smith #define __FUNC__ "SNESSetHessian" /* ADIC Ignore */ 10199b94acceSBarry Smith /*@C 10209b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 1021044dda88SLois Curfman McInnes location to store the matrix. 10229b94acceSBarry Smith 10239b94acceSBarry Smith Input Parameters: 10249b94acceSBarry Smith . snes - the SNES context 10259b94acceSBarry Smith . A - Hessian matrix 10269b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 10279b94acceSBarry Smith . func - Jacobian evaluation routine 1028313e4042SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 1029313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 10309b94acceSBarry Smith 10319b94acceSBarry Smith Calling sequence of func: 1032313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 10339b94acceSBarry Smith 10349b94acceSBarry Smith . x - input vector 10359b94acceSBarry Smith . A - Hessian matrix 10369b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1037ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1038ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 10392cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Hessian context 10409b94acceSBarry Smith 10419b94acceSBarry Smith Notes: 1042dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10432cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1044ac21db08SLois Curfman McInnes 10459b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 10469b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 10479b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10489b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10499b94acceSBarry Smith throughout the global iterations. 10509b94acceSBarry Smith 10519b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 10529b94acceSBarry Smith 1053ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 10549b94acceSBarry Smith @*/ 10559b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 10569b94acceSBarry Smith MatStructure*,void*),void *ctx) 10579b94acceSBarry Smith { 105877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 105974679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1060e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 10619b94acceSBarry Smith snes->computejacobian = func; 10629b94acceSBarry Smith snes->jacP = ctx; 10639b94acceSBarry Smith snes->jacobian = A; 10649b94acceSBarry Smith snes->jacobian_pre = B; 10659b94acceSBarry Smith return 0; 10669b94acceSBarry Smith } 10679b94acceSBarry Smith 10685615d1e5SSatish Balay #undef __FUNC__ 1069*5eea60f9SBarry Smith #define __FUNC__ "SNESGetHessian" /* ADIC Ignore */ 107062fef451SLois Curfman McInnes /*@ 107162fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 107262fef451SLois Curfman McInnes provided context for evaluating the Hessian. 107362fef451SLois Curfman McInnes 107462fef451SLois Curfman McInnes Input Parameter: 107562fef451SLois Curfman McInnes . snes - the nonlinear solver context 107662fef451SLois Curfman McInnes 107762fef451SLois Curfman McInnes Output Parameters: 107862fef451SLois Curfman McInnes . A - location to stash Hessian matrix (or PETSC_NULL) 107962fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 108062fef451SLois Curfman McInnes . ctx - location to stash Hessian ctx (or PETSC_NULL) 108162fef451SLois Curfman McInnes 108262fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 108362fef451SLois Curfman McInnes @*/ 108462fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 108562fef451SLois Curfman McInnes { 108674679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1087e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 108862fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 108962fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 109062fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 109162fef451SLois Curfman McInnes return 0; 109262fef451SLois Curfman McInnes } 109362fef451SLois Curfman McInnes 10949b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 10959b94acceSBarry Smith 10965615d1e5SSatish Balay #undef __FUNC__ 10975615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 10989b94acceSBarry Smith /*@ 10999b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1100272ac6f2SLois Curfman McInnes of a nonlinear solver. 11019b94acceSBarry Smith 11029b94acceSBarry Smith Input Parameter: 11039b94acceSBarry Smith . snes - the SNES context 11048ddd3da0SLois Curfman McInnes . x - the solution vector 11059b94acceSBarry Smith 1106272ac6f2SLois Curfman McInnes Notes: 1107272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1108272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1109272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1110272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1111272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1112272ac6f2SLois Curfman McInnes 11139b94acceSBarry Smith .keywords: SNES, nonlinear, setup 11149b94acceSBarry Smith 11159b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 11169b94acceSBarry Smith @*/ 11178ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 11189b94acceSBarry Smith { 1119272ac6f2SLois Curfman McInnes int ierr, flg; 112077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 112177c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 11228ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 11239b94acceSBarry Smith 1124c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1125c1f60f51SBarry Smith /* 1126c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1127dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1128c1f60f51SBarry Smith */ 1129c1f60f51SBarry Smith if (flg) { 1130c1f60f51SBarry Smith Mat J; 1131c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1132c1f60f51SBarry Smith PLogObjectParent(snes,J); 1133c1f60f51SBarry Smith snes->mfshell = J; 1134c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1135c1f60f51SBarry Smith snes->jacobian = J; 1136c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1137c1f60f51SBarry Smith } 1138c1f60f51SBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1139c1f60f51SBarry Smith snes->jacobian = J; 1140c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1141e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free operator option"); 1142c1f60f51SBarry Smith } 1143272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1144c1f60f51SBarry Smith /* 1145dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1146c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1147c1f60f51SBarry Smith */ 114831615d04SLois Curfman McInnes if (flg) { 1149272ac6f2SLois Curfman McInnes Mat J; 1150272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1151272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1152272ac6f2SLois Curfman McInnes snes->mfshell = J; 115383e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 115483e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 115594a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 115683e56afdSLois Curfman McInnes } 115783e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 115883e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 115994a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1160e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free option"); 1161272ac6f2SLois Curfman McInnes } 11629b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1163e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1164e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1165e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first"); 1166e3372554SBarry Smith if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector"); 1167a703fe33SLois Curfman McInnes 1168a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 116940191667SLois Curfman McInnes if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1170a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1171a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1172a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1173a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1174a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1175a703fe33SLois Curfman McInnes } 11769b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1177e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first"); 1178e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first"); 117937fcc0dbSBarry Smith if (!snes->computeumfunction) 1180e3372554SBarry Smith SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first"); 1181e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first"); 1182e3372554SBarry Smith } else SETERRQ(1,0,"Unknown method class"); 1183a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1184a703fe33SLois Curfman McInnes snes->setup_called = 1; 1185a703fe33SLois Curfman McInnes return 0; 11869b94acceSBarry Smith } 11879b94acceSBarry Smith 11885615d1e5SSatish Balay #undef __FUNC__ 1189*5eea60f9SBarry Smith #define __FUNC__ "SNESDestroy" /* ADIC Ignore */ 11909b94acceSBarry Smith /*@C 11919b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 11929b94acceSBarry Smith with SNESCreate(). 11939b94acceSBarry Smith 11949b94acceSBarry Smith Input Parameter: 11959b94acceSBarry Smith . snes - the SNES context 11969b94acceSBarry Smith 11979b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 11989b94acceSBarry Smith 119963a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 12009b94acceSBarry Smith @*/ 12019b94acceSBarry Smith int SNESDestroy(SNES snes) 12029b94acceSBarry Smith { 12039b94acceSBarry Smith int ierr; 120477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12059b94acceSBarry Smith ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr); 12060452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 12079b94acceSBarry Smith if (snes->mfshell) MatDestroy(snes->mfshell); 12089b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 12093e2e452bSBarry Smith if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 12106f24a144SLois Curfman McInnes if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 12119b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 12120452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 12139b94acceSBarry Smith return 0; 12149b94acceSBarry Smith } 12159b94acceSBarry Smith 12169b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 12179b94acceSBarry Smith 12185615d1e5SSatish Balay #undef __FUNC__ 12195615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 12209b94acceSBarry Smith /*@ 1221d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 12229b94acceSBarry Smith 12239b94acceSBarry Smith Input Parameters: 12249b94acceSBarry Smith . snes - the SNES context 122533174efeSLois Curfman McInnes . atol - absolute convergence tolerance 122633174efeSLois Curfman McInnes . rtol - relative convergence tolerance 122733174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 122833174efeSLois Curfman McInnes of the change in the solution between steps 122933174efeSLois Curfman McInnes . maxit - maximum number of iterations 123033174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 12319b94acceSBarry Smith 123233174efeSLois Curfman McInnes Options Database Keys: 123333174efeSLois Curfman McInnes $ -snes_atol <atol> 123433174efeSLois Curfman McInnes $ -snes_rtol <rtol> 123533174efeSLois Curfman McInnes $ -snes_stol <stol> 123633174efeSLois Curfman McInnes $ -snes_max_it <maxit> 123733174efeSLois Curfman McInnes $ -snes_max_funcs <maxf> 12389b94acceSBarry Smith 1239d7a720efSLois Curfman McInnes Notes: 12409b94acceSBarry Smith The default maximum number of iterations is 50. 12419b94acceSBarry Smith The default maximum number of function evaluations is 1000. 12429b94acceSBarry Smith 124333174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 12449b94acceSBarry Smith 1245d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 12469b94acceSBarry Smith @*/ 1247d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 12489b94acceSBarry Smith { 124977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1250d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1251d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1252d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1253d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1254d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 12559b94acceSBarry Smith return 0; 12569b94acceSBarry Smith } 12579b94acceSBarry Smith 12585615d1e5SSatish Balay #undef __FUNC__ 12595615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 12609b94acceSBarry Smith /*@ 126133174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 126233174efeSLois Curfman McInnes 126333174efeSLois Curfman McInnes Input Parameters: 126433174efeSLois Curfman McInnes . snes - the SNES context 126533174efeSLois Curfman McInnes . atol - absolute convergence tolerance 126633174efeSLois Curfman McInnes . rtol - relative convergence tolerance 126733174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 126833174efeSLois Curfman McInnes of the change in the solution between steps 126933174efeSLois Curfman McInnes . maxit - maximum number of iterations 127033174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 127133174efeSLois Curfman McInnes 127233174efeSLois Curfman McInnes Notes: 127333174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 127433174efeSLois Curfman McInnes 127533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 127633174efeSLois Curfman McInnes 127733174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 127833174efeSLois Curfman McInnes @*/ 127933174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 128033174efeSLois Curfman McInnes { 128133174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 128233174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 128333174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 128433174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 128533174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 128633174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 128733174efeSLois Curfman McInnes return 0; 128833174efeSLois Curfman McInnes } 128933174efeSLois Curfman McInnes 12905615d1e5SSatish Balay #undef __FUNC__ 12915615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 129233174efeSLois Curfman McInnes /*@ 12939b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 12949b94acceSBarry Smith 12959b94acceSBarry Smith Input Parameters: 12969b94acceSBarry Smith . snes - the SNES context 12979b94acceSBarry Smith . tol - tolerance 12989b94acceSBarry Smith 12999b94acceSBarry Smith Options Database Key: 1300d7a720efSLois Curfman McInnes $ -snes_trtol <tol> 13019b94acceSBarry Smith 13029b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 13039b94acceSBarry Smith 1304d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 13059b94acceSBarry Smith @*/ 13069b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 13079b94acceSBarry Smith { 130877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13099b94acceSBarry Smith snes->deltatol = tol; 13109b94acceSBarry Smith return 0; 13119b94acceSBarry Smith } 13129b94acceSBarry Smith 13135615d1e5SSatish Balay #undef __FUNC__ 13145615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 13159b94acceSBarry Smith /*@ 131677c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 13179b94acceSBarry Smith for unconstrained minimization solvers. 13189b94acceSBarry Smith 13199b94acceSBarry Smith Input Parameters: 13209b94acceSBarry Smith . snes - the SNES context 13219b94acceSBarry Smith . ftol - minimum function tolerance 13229b94acceSBarry Smith 13239b94acceSBarry Smith Options Database Key: 1324d7a720efSLois Curfman McInnes $ -snes_fmin <ftol> 13259b94acceSBarry Smith 13269b94acceSBarry Smith Note: 132777c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 13289b94acceSBarry Smith methods only. 13299b94acceSBarry Smith 13309b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 13319b94acceSBarry Smith 1332d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 13339b94acceSBarry Smith @*/ 133477c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 13359b94acceSBarry Smith { 133677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13379b94acceSBarry Smith snes->fmin = ftol; 13389b94acceSBarry Smith return 0; 13399b94acceSBarry Smith } 13409b94acceSBarry Smith 13419b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 13429b94acceSBarry Smith 13435615d1e5SSatish Balay #undef __FUNC__ 1344*5eea60f9SBarry Smith #define __FUNC__ "SNESSetMonitor" /* ADIC Ignore */ 13459b94acceSBarry Smith /*@C 13469b94acceSBarry Smith SNESSetMonitor - Sets the function that is to be used at every 13479b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 13489b94acceSBarry Smith progress. 13499b94acceSBarry Smith 13509b94acceSBarry Smith Input Parameters: 13519b94acceSBarry Smith . snes - the SNES context 13529b94acceSBarry Smith . func - monitoring routine 1353044dda88SLois Curfman McInnes . mctx - [optional] user-defined context for private data for the 1354044dda88SLois Curfman McInnes monitor routine (may be PETSC_NULL) 13559b94acceSBarry Smith 13569b94acceSBarry Smith Calling sequence of func: 1357682d7d0cSBarry Smith int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 13589b94acceSBarry Smith 13599b94acceSBarry Smith $ snes - the SNES context 13609b94acceSBarry Smith $ its - iteration number 1361044dda88SLois Curfman McInnes $ mctx - [optional] monitoring context 13629b94acceSBarry Smith $ 13639b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13649b94acceSBarry Smith $ norm - 2-norm function value (may be estimated) 13659b94acceSBarry Smith $ 13669b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13679b94acceSBarry Smith $ norm - 2-norm gradient value (may be estimated) 13689b94acceSBarry Smith 13699665c990SLois Curfman McInnes Options Database Keys: 13709665c990SLois Curfman McInnes $ -snes_monitor : sets SNESDefaultMonitor() 13719665c990SLois Curfman McInnes $ -snes_xmonitor : sets line graph monitor, 13729665c990SLois Curfman McInnes $ uses SNESLGMonitorCreate() 13739665c990SLois Curfman McInnes $ -snes_cancelmonitors : cancels all monitors that have 13749665c990SLois Curfman McInnes $ been hardwired into a code by 13759665c990SLois Curfman McInnes $ calls to SNESSetMonitor(), but 13769665c990SLois Curfman McInnes $ does not cancel those set via 13779665c990SLois Curfman McInnes $ the options database. 13789665c990SLois Curfman McInnes 13799665c990SLois Curfman McInnes 1380639f9d9dSBarry Smith Notes: 13816bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 13826bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 13836bc08f3fSLois Curfman McInnes order in which they were set. 1384639f9d9dSBarry Smith 13859b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 13869b94acceSBarry Smith 13879b94acceSBarry Smith .seealso: SNESDefaultMonitor() 13889b94acceSBarry Smith @*/ 138974679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 13909b94acceSBarry Smith { 1391639f9d9dSBarry Smith if (!func) { 1392639f9d9dSBarry Smith snes->numbermonitors = 0; 1393639f9d9dSBarry Smith return 0; 1394639f9d9dSBarry Smith } 1395639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1396e3372554SBarry Smith SETERRQ(1,0,"Too many monitors set"); 1397639f9d9dSBarry Smith } 1398639f9d9dSBarry Smith 1399639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1400639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 14019b94acceSBarry Smith return 0; 14029b94acceSBarry Smith } 14039b94acceSBarry Smith 14045615d1e5SSatish Balay #undef __FUNC__ 1405*5eea60f9SBarry Smith #define __FUNC__ "SNESSetConvergenceTest" /* ADIC Ignore */ 14069b94acceSBarry Smith /*@C 14079b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 14089b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 14099b94acceSBarry Smith 14109b94acceSBarry Smith Input Parameters: 14119b94acceSBarry Smith . snes - the SNES context 14129b94acceSBarry Smith . func - routine to test for convergence 1413044dda88SLois Curfman McInnes . cctx - [optional] context for private data for the convergence routine 1414044dda88SLois Curfman McInnes (may be PETSC_NULL) 14159b94acceSBarry Smith 14169b94acceSBarry Smith Calling sequence of func: 14179b94acceSBarry Smith int func (SNES snes,double xnorm,double gnorm, 14189b94acceSBarry Smith double f,void *cctx) 14199b94acceSBarry Smith 14209b94acceSBarry Smith $ snes - the SNES context 1421044dda88SLois Curfman McInnes $ cctx - [optional] convergence context 14229b94acceSBarry Smith $ xnorm - 2-norm of current iterate 14239b94acceSBarry Smith $ 14249b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 14259b94acceSBarry Smith $ gnorm - 2-norm of current step 14269b94acceSBarry Smith $ f - 2-norm of function 14279b94acceSBarry Smith $ 14289b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 14299b94acceSBarry Smith $ gnorm - 2-norm of current gradient 14309b94acceSBarry Smith $ f - function value 14319b94acceSBarry Smith 14329b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 14339b94acceSBarry Smith 143440191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 143540191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 14369b94acceSBarry Smith @*/ 143774679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 14389b94acceSBarry Smith { 14399b94acceSBarry Smith (snes)->converged = func; 14409b94acceSBarry Smith (snes)->cnvP = cctx; 14419b94acceSBarry Smith return 0; 14429b94acceSBarry Smith } 14439b94acceSBarry Smith 14445615d1e5SSatish Balay #undef __FUNC__ 1445*5eea60f9SBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" /* ADIC Ignore */ 1446c9005455SLois Curfman McInnes /*@ 1447c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1448c9005455SLois Curfman McInnes 1449c9005455SLois Curfman McInnes Input Parameters: 1450c9005455SLois Curfman McInnes . snes - iterative context obtained from SNESCreate() 1451c9005455SLois Curfman McInnes . a - array to hold history 1452c9005455SLois Curfman McInnes . na - size of a 1453c9005455SLois Curfman McInnes 1454c9005455SLois Curfman McInnes Notes: 1455c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1456c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1457c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1458c9005455SLois Curfman McInnes at each step. 1459c9005455SLois Curfman McInnes 1460c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1461c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1462c9005455SLois Curfman McInnes during the section of code that is being timed. 1463c9005455SLois Curfman McInnes 1464c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1465c9005455SLois Curfman McInnes @*/ 1466c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na) 1467c9005455SLois Curfman McInnes { 1468c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1469c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1470c9005455SLois Curfman McInnes snes->conv_hist = a; 1471c9005455SLois Curfman McInnes snes->conv_hist_size = na; 1472c9005455SLois Curfman McInnes return 0; 1473c9005455SLois Curfman McInnes } 1474c9005455SLois Curfman McInnes 1475c9005455SLois Curfman McInnes #undef __FUNC__ 14765615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 14779b94acceSBarry Smith /* 14789b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 14799b94acceSBarry Smith positive parameter delta. 14809b94acceSBarry Smith 14819b94acceSBarry Smith Input Parameters: 14829b94acceSBarry Smith . snes - the SNES context 14839b94acceSBarry Smith . y - approximate solution of linear system 14849b94acceSBarry Smith . fnorm - 2-norm of current function 14859b94acceSBarry Smith . delta - trust region size 14869b94acceSBarry Smith 14879b94acceSBarry Smith Output Parameters: 14889b94acceSBarry Smith . gpnorm - predicted function norm at the new point, assuming local 14899b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 14909b94acceSBarry Smith region, and exceeds zero otherwise. 14919b94acceSBarry Smith . ynorm - 2-norm of the step 14929b94acceSBarry Smith 14939b94acceSBarry Smith Note: 149440191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 14959b94acceSBarry Smith is set to be the maximum allowable step size. 14969b94acceSBarry Smith 14979b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 14989b94acceSBarry Smith */ 14999b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 15009b94acceSBarry Smith double *gpnorm,double *ynorm) 15019b94acceSBarry Smith { 15029b94acceSBarry Smith double norm; 15039b94acceSBarry Smith Scalar cnorm; 1504cddf8d76SBarry Smith VecNorm(y,NORM_2, &norm ); 15059b94acceSBarry Smith if (norm > *delta) { 15069b94acceSBarry Smith norm = *delta/norm; 15079b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 15089b94acceSBarry Smith cnorm = norm; 15099b94acceSBarry Smith VecScale( &cnorm, y ); 15109b94acceSBarry Smith *ynorm = *delta; 15119b94acceSBarry Smith } else { 15129b94acceSBarry Smith *gpnorm = 0.0; 15139b94acceSBarry Smith *ynorm = norm; 15149b94acceSBarry Smith } 15159b94acceSBarry Smith return 0; 15169b94acceSBarry Smith } 15179b94acceSBarry Smith 15185615d1e5SSatish Balay #undef __FUNC__ 15195615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 15209b94acceSBarry Smith /*@ 15219b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 152263a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 15239b94acceSBarry Smith 15249b94acceSBarry Smith Input Parameter: 15259b94acceSBarry Smith . snes - the SNES context 15268ddd3da0SLois Curfman McInnes . x - the solution vector 15279b94acceSBarry Smith 15289b94acceSBarry Smith Output Parameter: 15299b94acceSBarry Smith its - number of iterations until termination 15309b94acceSBarry Smith 15318ddd3da0SLois Curfman McInnes Note: 15328ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 15338ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 15348ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 15358ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 15368ddd3da0SLois Curfman McInnes 15379b94acceSBarry Smith .keywords: SNES, nonlinear, solve 15389b94acceSBarry Smith 153963a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 15409b94acceSBarry Smith @*/ 15418ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 15429b94acceSBarry Smith { 15433c7409f5SSatish Balay int ierr, flg; 1544052efed2SBarry Smith 154577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 154674679c65SBarry Smith PetscValidIntPointer(its); 1547c4fc05e7SBarry Smith if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1548c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 15499b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 1550c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 15519b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 15529b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 15533c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 15546d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 15559b94acceSBarry Smith return 0; 15569b94acceSBarry Smith } 15579b94acceSBarry Smith 15589b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 155937fcc0dbSBarry Smith static NRList *__SNESList = 0; 15609b94acceSBarry Smith 15615615d1e5SSatish Balay #undef __FUNC__ 15625615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 15639b94acceSBarry Smith /*@ 15644b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 15659b94acceSBarry Smith 15669b94acceSBarry Smith Input Parameters: 15679b94acceSBarry Smith . snes - the SNES context 15689b94acceSBarry Smith . method - a known method 15699b94acceSBarry Smith 1570ae12b187SLois Curfman McInnes Options Database Command: 1571ae12b187SLois Curfman McInnes $ -snes_type <method> 1572ae12b187SLois Curfman McInnes $ Use -help for a list of available methods 1573ae12b187SLois Curfman McInnes $ (for instance, ls or tr) 1574ae12b187SLois Curfman McInnes 15759b94acceSBarry Smith Notes: 15769b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 15779b94acceSBarry Smith $ Systems of nonlinear equations: 157840191667SLois Curfman McInnes $ SNES_EQ_LS - Newton's method with line search 157940191667SLois Curfman McInnes $ SNES_EQ_TR - Newton's method with trust region 15809b94acceSBarry Smith $ Unconstrained minimization: 158140191667SLois Curfman McInnes $ SNES_UM_TR - Newton's method with trust region 158240191667SLois Curfman McInnes $ SNES_UM_LS - Newton's method with line search 15839b94acceSBarry Smith 1584ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1585ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1586ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1587ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1588ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1589ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1590ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1591ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1592ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1593ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1594a703fe33SLois Curfman McInnes 1595f536c699SSatish Balay .keywords: SNES, set, method 15969b94acceSBarry Smith @*/ 15974b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 15989b94acceSBarry Smith { 159984cb2905SBarry Smith int ierr; 16009b94acceSBarry Smith int (*r)(SNES); 160177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16027d1a2b2bSBarry Smith if (snes->type == (int) method) return 0; 16037d1a2b2bSBarry Smith 16049b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 160584cb2905SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1606e3372554SBarry Smith if (!__SNESList) {SETERRQ(1,0,"Could not get methods");} 160737fcc0dbSBarry Smith r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1608e3372554SBarry Smith if (!r) {SETERRQ(1,0,"Unknown method");} 16090452661fSBarry Smith if (snes->data) PetscFree(snes->data); 16109b94acceSBarry Smith snes->set_method_called = 1; 16119b94acceSBarry Smith return (*r)(snes); 16129b94acceSBarry Smith } 16139b94acceSBarry Smith 16149b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16155615d1e5SSatish Balay #undef __FUNC__ 1616*5eea60f9SBarry Smith #define __FUNC__ "SNESRegister" /* ADIC Ignore */ 16179b94acceSBarry Smith /*@C 16189b94acceSBarry Smith SNESRegister - Adds the method to the nonlinear solver package, given 16194b0e389bSBarry Smith a function pointer and a nonlinear solver name of the type SNESType. 16209b94acceSBarry Smith 16219b94acceSBarry Smith Input Parameters: 16222d872ea7SLois Curfman McInnes . name - either a predefined name such as SNES_EQ_LS, or SNES_NEW 16232d872ea7SLois Curfman McInnes to indicate a new user-defined solver 162440191667SLois Curfman McInnes . sname - corresponding string for name 16259b94acceSBarry Smith . create - routine to create method context 16269b94acceSBarry Smith 162784cb2905SBarry Smith Output Parameter: 162884cb2905SBarry Smith . oname - type associated with this new solver 162984cb2905SBarry Smith 16302d872ea7SLois Curfman McInnes Notes: 16312d872ea7SLois Curfman McInnes Multiple user-defined nonlinear solvers can be added by calling 16322d872ea7SLois Curfman McInnes SNESRegister() with the input parameter "name" set to be SNES_NEW; 16332d872ea7SLois Curfman McInnes each call will return a unique solver type in the output 16342d872ea7SLois Curfman McInnes parameter "oname". 16352d872ea7SLois Curfman McInnes 16369b94acceSBarry Smith .keywords: SNES, nonlinear, register 16379b94acceSBarry Smith 16389b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy() 16399b94acceSBarry Smith @*/ 164084cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES)) 16419b94acceSBarry Smith { 16429b94acceSBarry Smith int ierr; 164384cb2905SBarry Smith static int numberregistered = 0; 164484cb2905SBarry Smith 1645d252947aSBarry Smith if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++); 164684cb2905SBarry Smith 164784cb2905SBarry Smith if (oname) *oname = name; 164837fcc0dbSBarry Smith if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 164984cb2905SBarry Smith NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create ); 16509b94acceSBarry Smith return 0; 16519b94acceSBarry Smith } 1652a847f771SSatish Balay 16539b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16545615d1e5SSatish Balay #undef __FUNC__ 1655*5eea60f9SBarry Smith #define __FUNC__ "SNESRegisterDestroy" /* ADIC Ignore */ 16569b94acceSBarry Smith /*@C 16579b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 16589b94acceSBarry Smith registered by SNESRegister(). 16599b94acceSBarry Smith 16609b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 16619b94acceSBarry Smith 16629b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 16639b94acceSBarry Smith @*/ 16649b94acceSBarry Smith int SNESRegisterDestroy() 16659b94acceSBarry Smith { 166637fcc0dbSBarry Smith if (__SNESList) { 166737fcc0dbSBarry Smith NRDestroy( __SNESList ); 166837fcc0dbSBarry Smith __SNESList = 0; 16699b94acceSBarry Smith } 167084cb2905SBarry Smith SNESRegisterAllCalled = 0; 16719b94acceSBarry Smith return 0; 16729b94acceSBarry Smith } 16739b94acceSBarry Smith 16745615d1e5SSatish Balay #undef __FUNC__ 1675*5eea60f9SBarry Smith #define __FUNC__ "SNESGetTypeFromOptions_Private" /* ADIC Ignore */ 16769b94acceSBarry Smith /* 16774b0e389bSBarry Smith SNESGetTypeFromOptions_Private - Sets the selected method from the 16789b94acceSBarry Smith options database. 16799b94acceSBarry Smith 16809b94acceSBarry Smith Input Parameter: 16819b94acceSBarry Smith . ctx - the SNES context 16829b94acceSBarry Smith 16839b94acceSBarry Smith Output Parameter: 16849b94acceSBarry Smith . method - solver method 16859b94acceSBarry Smith 16869b94acceSBarry Smith Returns: 16879b94acceSBarry Smith Returns 1 if the method is found; 0 otherwise. 16889b94acceSBarry Smith 16899b94acceSBarry Smith Options Database Key: 16904b0e389bSBarry Smith $ -snes_type method 16919b94acceSBarry Smith */ 1692052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 16939b94acceSBarry Smith { 1694052efed2SBarry Smith int ierr; 16959b94acceSBarry Smith char sbuf[50]; 16965baf8537SBarry Smith 1697052efed2SBarry Smith ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1698052efed2SBarry Smith if (*flg) { 1699052efed2SBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 170037fcc0dbSBarry Smith *method = (SNESType)NRFindID( __SNESList, sbuf ); 17019b94acceSBarry Smith } 17029b94acceSBarry Smith return 0; 17039b94acceSBarry Smith } 17049b94acceSBarry Smith 17055615d1e5SSatish Balay #undef __FUNC__ 1706*5eea60f9SBarry Smith #define __FUNC__ "SNESGetType" /* ADIC Ignore */ 17079b94acceSBarry Smith /*@C 17089a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 17099b94acceSBarry Smith 17109b94acceSBarry Smith Input Parameter: 17114b0e389bSBarry Smith . snes - nonlinear solver context 17129b94acceSBarry Smith 17139b94acceSBarry Smith Output Parameter: 17149a28b0a6SLois Curfman McInnes . method - SNES method (or use PETSC_NULL) 17159a28b0a6SLois Curfman McInnes . name - name of SNES method (or use PETSC_NULL) 17169b94acceSBarry Smith 17179b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 17189b94acceSBarry Smith @*/ 17194b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name) 17209b94acceSBarry Smith { 172106a9b0adSLois Curfman McInnes int ierr; 172237fcc0dbSBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 17234b0e389bSBarry Smith if (method) *method = (SNESType) snes->type; 172437fcc0dbSBarry Smith if (name) *name = NRFindName( __SNESList, (int) snes->type ); 17259b94acceSBarry Smith return 0; 17269b94acceSBarry Smith } 17279b94acceSBarry Smith 17289b94acceSBarry Smith #include <stdio.h> 17295615d1e5SSatish Balay #undef __FUNC__ 1730*5eea60f9SBarry Smith #define __FUNC__ "SNESPrintTypes_Private" /* ADIC Ignore */ 17319b94acceSBarry Smith /* 17324b0e389bSBarry Smith SNESPrintTypes_Private - Prints the SNES methods available from the 17339b94acceSBarry Smith options database. 17349b94acceSBarry Smith 17359b94acceSBarry Smith Input Parameters: 173633455573SSatish Balay . comm - communicator (usually MPI_COMM_WORLD) 17379b94acceSBarry Smith . prefix - prefix (usually "-") 17384b0e389bSBarry Smith . name - the options database name (by default "snes_type") 17399b94acceSBarry Smith */ 174033455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 17419b94acceSBarry Smith { 17429b94acceSBarry Smith FuncList *entry; 174337fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 174437fcc0dbSBarry Smith entry = __SNESList->head; 174577c4ece6SBarry Smith PetscPrintf(comm," %s%s (one of)",prefix,name); 17469b94acceSBarry Smith while (entry) { 174777c4ece6SBarry Smith PetscPrintf(comm," %s",entry->name); 17489b94acceSBarry Smith entry = entry->next; 17499b94acceSBarry Smith } 175077c4ece6SBarry Smith PetscPrintf(comm,"\n"); 17519b94acceSBarry Smith return 0; 17529b94acceSBarry Smith } 17539b94acceSBarry Smith 17545615d1e5SSatish Balay #undef __FUNC__ 1755*5eea60f9SBarry Smith #define __FUNC__ "SNESGetSolution" /* ADIC Ignore */ 17569b94acceSBarry Smith /*@C 17579b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 17589b94acceSBarry Smith stored. 17599b94acceSBarry Smith 17609b94acceSBarry Smith Input Parameter: 17619b94acceSBarry Smith . snes - the SNES context 17629b94acceSBarry Smith 17639b94acceSBarry Smith Output Parameter: 17649b94acceSBarry Smith . x - the solution 17659b94acceSBarry Smith 17669b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 17679b94acceSBarry Smith 1768a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 17699b94acceSBarry Smith @*/ 17709b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 17719b94acceSBarry Smith { 177277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17739b94acceSBarry Smith *x = snes->vec_sol_always; 17749b94acceSBarry Smith return 0; 17759b94acceSBarry Smith } 17769b94acceSBarry Smith 17775615d1e5SSatish Balay #undef __FUNC__ 1778*5eea60f9SBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" /* ADIC Ignore */ 17799b94acceSBarry Smith /*@C 17809b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 17819b94acceSBarry Smith stored. 17829b94acceSBarry Smith 17839b94acceSBarry Smith Input Parameter: 17849b94acceSBarry Smith . snes - the SNES context 17859b94acceSBarry Smith 17869b94acceSBarry Smith Output Parameter: 17879b94acceSBarry Smith . x - the solution update 17889b94acceSBarry Smith 17899b94acceSBarry Smith Notes: 17909b94acceSBarry Smith This vector is implementation dependent. 17919b94acceSBarry Smith 17929b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 17939b94acceSBarry Smith 17949b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 17959b94acceSBarry Smith @*/ 17969b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 17979b94acceSBarry Smith { 179877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17999b94acceSBarry Smith *x = snes->vec_sol_update_always; 18009b94acceSBarry Smith return 0; 18019b94acceSBarry Smith } 18029b94acceSBarry Smith 18035615d1e5SSatish Balay #undef __FUNC__ 1804*5eea60f9SBarry Smith #define __FUNC__ "SNESGetFunction" /* ADIC Ignore */ 18059b94acceSBarry Smith /*@C 18063638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 18079b94acceSBarry Smith 18089b94acceSBarry Smith Input Parameter: 18099b94acceSBarry Smith . snes - the SNES context 18109b94acceSBarry Smith 18119b94acceSBarry Smith Output Parameter: 18123638b69dSLois Curfman McInnes . r - the function 18139b94acceSBarry Smith 18149b94acceSBarry Smith Notes: 18159b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 18169b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 18179b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 18189b94acceSBarry Smith 1819a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 18209b94acceSBarry Smith 182131615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 182231615d04SLois Curfman McInnes SNESGetGradient() 18239b94acceSBarry Smith @*/ 18249b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 18259b94acceSBarry Smith { 182677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1827e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0, 18287c792091SSatish Balay "For SNES_NONLINEAR_EQUATIONS only"); 18299b94acceSBarry Smith *r = snes->vec_func_always; 18309b94acceSBarry Smith return 0; 18319b94acceSBarry Smith } 18329b94acceSBarry Smith 18335615d1e5SSatish Balay #undef __FUNC__ 1834*5eea60f9SBarry Smith #define __FUNC__ "SNESGetGradient" /* ADIC Ignore */ 18359b94acceSBarry Smith /*@C 18363638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 18379b94acceSBarry Smith 18389b94acceSBarry Smith Input Parameter: 18399b94acceSBarry Smith . snes - the SNES context 18409b94acceSBarry Smith 18419b94acceSBarry Smith Output Parameter: 18429b94acceSBarry Smith . r - the gradient 18439b94acceSBarry Smith 18449b94acceSBarry Smith Notes: 18459b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 18469b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18479b94acceSBarry Smith SNESGetFunction(). 18489b94acceSBarry Smith 18499b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 18509b94acceSBarry Smith 185131615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 18529b94acceSBarry Smith @*/ 18539b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 18549b94acceSBarry Smith { 185577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1856e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 185763c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 18589b94acceSBarry Smith *r = snes->vec_func_always; 18599b94acceSBarry Smith return 0; 18609b94acceSBarry Smith } 18619b94acceSBarry Smith 18625615d1e5SSatish Balay #undef __FUNC__ 1863*5eea60f9SBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" /* ADIC Ignore */ 18649b94acceSBarry Smith /*@ 18659b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 18669b94acceSBarry Smith unconstrained minimization problems. 18679b94acceSBarry Smith 18689b94acceSBarry Smith Input Parameter: 18699b94acceSBarry Smith . snes - the SNES context 18709b94acceSBarry Smith 18719b94acceSBarry Smith Output Parameter: 18729b94acceSBarry Smith . r - the function 18739b94acceSBarry Smith 18749b94acceSBarry Smith Notes: 18759b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 18769b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18779b94acceSBarry Smith SNESGetFunction(). 18789b94acceSBarry Smith 18799b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 18809b94acceSBarry Smith 188131615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 18829b94acceSBarry Smith @*/ 18839b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 18849b94acceSBarry Smith { 188577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 188674679c65SBarry Smith PetscValidScalarPointer(r); 1887e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 188863c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 18899b94acceSBarry Smith *r = snes->fc; 18909b94acceSBarry Smith return 0; 18919b94acceSBarry Smith } 18929b94acceSBarry Smith 18935615d1e5SSatish Balay #undef __FUNC__ 1894*5eea60f9SBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" /* ADIC Ignore */ 18953c7409f5SSatish Balay /*@C 18963c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1897639f9d9dSBarry Smith SNES options in the database. You must NOT include the - at the beginning of 1898639f9d9dSBarry Smith the prefix name. 18993c7409f5SSatish Balay 19003c7409f5SSatish Balay Input Parameter: 19013c7409f5SSatish Balay . snes - the SNES context 19023c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 19033c7409f5SSatish Balay 19043c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1905a86d99e1SLois Curfman McInnes 1906a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19073c7409f5SSatish Balay @*/ 19083c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 19093c7409f5SSatish Balay { 19103c7409f5SSatish Balay int ierr; 19113c7409f5SSatish Balay 191277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1913639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19143c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19153c7409f5SSatish Balay return 0; 19163c7409f5SSatish Balay } 19173c7409f5SSatish Balay 19185615d1e5SSatish Balay #undef __FUNC__ 1919*5eea60f9SBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" /* ADIC Ignore */ 19203c7409f5SSatish Balay /*@C 1921f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1922639f9d9dSBarry Smith SNES options in the database. You must NOT include the - at the beginning of 1923639f9d9dSBarry Smith the prefix name. 19243c7409f5SSatish Balay 19253c7409f5SSatish Balay Input Parameter: 19263c7409f5SSatish Balay . snes - the SNES context 19273c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 19283c7409f5SSatish Balay 19293c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1930a86d99e1SLois Curfman McInnes 1931a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 19323c7409f5SSatish Balay @*/ 19333c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 19343c7409f5SSatish Balay { 19353c7409f5SSatish Balay int ierr; 19363c7409f5SSatish Balay 193777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1938639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19393c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19403c7409f5SSatish Balay return 0; 19413c7409f5SSatish Balay } 19423c7409f5SSatish Balay 19435615d1e5SSatish Balay #undef __FUNC__ 1944*5eea60f9SBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" /* ADIC Ignore */ 1945c83590e2SSatish Balay /*@ 19463c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 19473c7409f5SSatish Balay SNES options in the database. 19483c7409f5SSatish Balay 19493c7409f5SSatish Balay Input Parameter: 19503c7409f5SSatish Balay . snes - the SNES context 19513c7409f5SSatish Balay 19523c7409f5SSatish Balay Output Parameter: 19533c7409f5SSatish Balay . prefix - pointer to the prefix string used 19543c7409f5SSatish Balay 19553c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 1956a86d99e1SLois Curfman McInnes 1957a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 19583c7409f5SSatish Balay @*/ 19593c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 19603c7409f5SSatish Balay { 19613c7409f5SSatish Balay int ierr; 19623c7409f5SSatish Balay 196377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1964639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19653c7409f5SSatish Balay return 0; 19663c7409f5SSatish Balay } 19673c7409f5SSatish Balay 19683c7409f5SSatish Balay 19699b94acceSBarry Smith 19709b94acceSBarry Smith 19719b94acceSBarry Smith 1972