1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*11ca99fdSLois Curfman McInnes static char vcid[] = "$Id: snes.c,v 1.128 1997/08/22 15:17:50 bsmith Exp curfman $"; 39b94acceSBarry Smith #endif 49b94acceSBarry Smith 570f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 6f5eb4b81SSatish Balay #include "src/sys/nreg.h" 79b94acceSBarry Smith #include "pinclude/pviewer.h" 89b94acceSBarry Smith #include <math.h> 99b94acceSBarry Smith 10052efed2SBarry Smith extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*); 1133455573SSatish Balay extern int SNESPrintTypes_Private(MPI_Comm,char*,char*); 129b94acceSBarry Smith 1384cb2905SBarry Smith int SNESRegisterAllCalled = 0; 1484cb2905SBarry Smith 155615d1e5SSatish Balay #undef __FUNC__ 16d4bb536fSBarry Smith #define __FUNC__ "SNESView" 179b94acceSBarry Smith /*@ 189b94acceSBarry Smith SNESView - Prints the SNES data structure. 199b94acceSBarry Smith 209b94acceSBarry Smith Input Parameters: 219b94acceSBarry Smith . SNES - the SNES context 229b94acceSBarry Smith . viewer - visualization context 239b94acceSBarry Smith 249b94acceSBarry Smith Options Database Key: 259b94acceSBarry Smith $ -snes_view : calls SNESView() at end of SNESSolve() 269b94acceSBarry Smith 279b94acceSBarry Smith Notes: 289b94acceSBarry Smith The available visualization contexts include 296d4a8577SBarry Smith $ VIEWER_STDOUT_SELF - standard output (default) 306d4a8577SBarry Smith $ VIEWER_STDOUT_WORLD - synchronized standard 319b94acceSBarry Smith $ output where only the first processor opens 329b94acceSBarry Smith $ the file. All other processors send their 339b94acceSBarry Smith $ data to the first processor to print. 349b94acceSBarry Smith 359b94acceSBarry Smith The user can open alternative vistualization contexts with 36dbb450caSBarry Smith $ ViewerFileOpenASCII() - output to a specified file 379b94acceSBarry Smith 389b94acceSBarry Smith .keywords: SNES, view 399b94acceSBarry Smith 40dbb450caSBarry Smith .seealso: ViewerFileOpenASCII() 419b94acceSBarry Smith @*/ 429b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 439b94acceSBarry Smith { 449b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 459b94acceSBarry Smith FILE *fd; 469b94acceSBarry Smith int ierr; 479b94acceSBarry Smith SLES sles; 489b94acceSBarry Smith char *method; 4919bcc07fSBarry Smith ViewerType vtype; 509b94acceSBarry Smith 5174679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5274679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5374679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5474679c65SBarry Smith 5519bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 5619bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 5790ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 5877c4ece6SBarry Smith PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 594b0e389bSBarry Smith SNESGetType(snes,PETSC_NULL,&method); 6077c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," method: %s\n",method); 619b94acceSBarry Smith if (snes->view) (*snes->view)((PetscObject)snes,viewer); 6277c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 639b94acceSBarry Smith " maximum iterations=%d, maximum function evaluations=%d\n", 649b94acceSBarry Smith snes->max_its,snes->max_funcs); 6577c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 669b94acceSBarry Smith " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 679b94acceSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol); 687a00f4a9SLois Curfman McInnes PetscFPrintf(snes->comm,fd, 697a00f4a9SLois Curfman McInnes " total number of linear solver iterations=%d\n",snes->linear_its); 70ae3c334cSLois Curfman McInnes PetscFPrintf(snes->comm,fd, 7168632166SLois Curfman McInnes " total number of function evaluations=%d\n",snes->nfuncs); 729b94acceSBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 7377c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 749b94acceSBarry Smith if (snes->ksp_ewconv) { 759b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 769b94acceSBarry Smith if (kctx) { 7777c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 789b94acceSBarry Smith " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 799b94acceSBarry Smith kctx->version); 8077c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 819b94acceSBarry Smith " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 829b94acceSBarry Smith kctx->rtol_max,kctx->threshold); 8377c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 849b94acceSBarry Smith kctx->gamma,kctx->alpha,kctx->alpha2); 859b94acceSBarry Smith } 869b94acceSBarry Smith } 87c30f285eSLois Curfman McInnes } else if (vtype == STRING_VIEWER) { 88c30f285eSLois Curfman McInnes SNESGetType(snes,PETSC_NULL,&method); 89c30f285eSLois Curfman McInnes ViewerStringSPrintf(viewer," %-3.3s",method); 9019bcc07fSBarry Smith } 919b94acceSBarry Smith SNESGetSLES(snes,&sles); 929b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 939b94acceSBarry Smith return 0; 949b94acceSBarry Smith } 959b94acceSBarry Smith 96639f9d9dSBarry Smith /* 97639f9d9dSBarry Smith We retain a list of functions that also take SNES command 98639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 99639f9d9dSBarry Smith */ 100639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 101639f9d9dSBarry Smith static int numberofsetfromoptions; 102639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 103639f9d9dSBarry Smith 1045615d1e5SSatish Balay #undef __FUNC__ 105d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker" 106639f9d9dSBarry Smith /*@ 107639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 108639f9d9dSBarry Smith 109639f9d9dSBarry Smith Input Parameter: 110639f9d9dSBarry Smith . snescheck - function that checks for options 111639f9d9dSBarry Smith 112639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 113639f9d9dSBarry Smith @*/ 114639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 115639f9d9dSBarry Smith { 116639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 117e3372554SBarry Smith SETERRQ(1,0,"Too many options checkers, only 5 allowed"); 118639f9d9dSBarry Smith } 119639f9d9dSBarry Smith 120639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 121639f9d9dSBarry Smith return 0; 122639f9d9dSBarry Smith } 123639f9d9dSBarry Smith 1245615d1e5SSatish Balay #undef __FUNC__ 1255615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1269b94acceSBarry Smith /*@ 127682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1289b94acceSBarry Smith 1299b94acceSBarry Smith Input Parameter: 1309b94acceSBarry Smith . snes - the SNES context 1319b94acceSBarry Smith 132*11ca99fdSLois Curfman McInnes Notes: 133*11ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 134*11ca99fdSLois Curfman McInnes the users manual. 13583e2fdc7SBarry Smith 1369b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1379b94acceSBarry Smith 138a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1399b94acceSBarry Smith @*/ 1409b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1419b94acceSBarry Smith { 1424b0e389bSBarry Smith SNESType method; 1439b94acceSBarry Smith double tmp; 1449b94acceSBarry Smith SLES sles; 14581f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1463c7409f5SSatish Balay int version = PETSC_DEFAULT; 1479b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1489b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1499b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1509b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1519b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1529b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1539b94acceSBarry Smith 15481f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 15581f57ec7SBarry Smith 15681f57ec7SBarry Smith 15777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 158e3372554SBarry Smith if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp"); 159052efed2SBarry Smith ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 160052efed2SBarry Smith if (flg) { 161052efed2SBarry Smith ierr = SNESSetType(snes,method); CHKERRQ(ierr); 1629b94acceSBarry Smith } 1635334005bSBarry Smith else if (!snes->set_method_called) { 1645334005bSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 16540191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 1665334005bSBarry Smith } 1675334005bSBarry Smith else { 16840191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 1695334005bSBarry Smith } 1705334005bSBarry Smith } 1715334005bSBarry Smith 1723c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1733c7409f5SSatish Balay if (flg) { SNESPrintHelp(snes); } 1743c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 175d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);} 1763c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 177d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);} 1783c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 179d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); } 1803c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 1813c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 182d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 183d7a720efSLois Curfman McInnes if (flg) { SNESSetTrustRegionTolerance(snes,tmp); } 184d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 185d7a720efSLois Curfman McInnes if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); } 1863c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 1873c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 188b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 189b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 190b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 191b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 192b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 193b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 194b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 19593c39befSBarry Smith 19693c39befSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr); 19793c39befSBarry Smith if (flg) {snes->converged = 0;} 19893c39befSBarry Smith 1999b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 2009b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 2015bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 2025bbad29bSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);} 2033c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 204639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 2053c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 206d31a3109SLois Curfman McInnes if (flg) { 207639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 208d31a3109SLois Curfman McInnes } 20981f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2103c7409f5SSatish Balay if (flg) { 21117699dbbSLois Curfman McInnes int rank = 0; 212d7e8b826SBarry Smith DrawLG lg; 21317699dbbSLois Curfman McInnes MPI_Initialized(&rank); 21417699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 21517699dbbSLois Curfman McInnes if (!rank) { 21681f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2179b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 218f8c826e1SBarry Smith snes->xmonitor = lg; 2199b94acceSBarry Smith PLogObjectParent(snes,lg); 2209b94acceSBarry Smith } 2219b94acceSBarry Smith } 2223c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2233c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2249b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2259b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 22694a424c1SBarry Smith PLogInfo(snes, 22731615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 22831615d04SLois Curfman McInnes } 22931615d04SLois Curfman McInnes else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 23031615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 23131615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 23294a424c1SBarry Smith PLogInfo(snes, 23331615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2349b94acceSBarry Smith } 235639f9d9dSBarry Smith 236639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 237639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 238639f9d9dSBarry Smith } 239639f9d9dSBarry Smith 2409b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2419b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 2429b94acceSBarry Smith if (!snes->setfromoptions) return 0; 2439b94acceSBarry Smith return (*snes->setfromoptions)(snes); 2449b94acceSBarry Smith } 2459b94acceSBarry Smith 2465615d1e5SSatish Balay #undef __FUNC__ 247d4bb536fSBarry Smith #define __FUNC__ "SNESPrintHelp" 2489b94acceSBarry Smith /*@ 2499b94acceSBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2509b94acceSBarry Smith 2519b94acceSBarry Smith Input Parameter: 2529b94acceSBarry Smith . snes - the SNES context 2539b94acceSBarry Smith 254a703fe33SLois Curfman McInnes Options Database Keys: 255a703fe33SLois Curfman McInnes $ -help, -h 256a703fe33SLois Curfman McInnes 2579b94acceSBarry Smith .keywords: SNES, nonlinear, help 2589b94acceSBarry Smith 259682d7d0cSBarry Smith .seealso: SNESSetFromOptions() 2609b94acceSBarry Smith @*/ 2619b94acceSBarry Smith int SNESPrintHelp(SNES snes) 2629b94acceSBarry Smith { 2636daaf66cSBarry Smith char p[64]; 2646daaf66cSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2656daaf66cSBarry Smith 26677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2676daaf66cSBarry Smith 2686daaf66cSBarry Smith PetscStrcpy(p,"-"); 2696daaf66cSBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2706daaf66cSBarry Smith 2716daaf66cSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2726daaf66cSBarry Smith 273d31a3109SLois Curfman McInnes PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 27433455573SSatish Balay SNESPrintTypes_Private(snes->comm,p,"snes_type"); 27577c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 276b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 277b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 278b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 279b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 280b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 281b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 2825bbad29bSBarry Smith PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 2835bbad29bSBarry Smith PetscPrintf(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 28483e2fdc7SBarry Smith PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 28583e2fdc7SBarry Smith residual norm at each iteration.\n",p); 28683e2fdc7SBarry Smith PetscPrintf(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 28783e2fdc7SBarry Smith residual norm for small residual norms. This is useful to conceal\n\ 28883e2fdc7SBarry Smith meaningless digits that may be different on different machines.\n",p); 289d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 290b18e04deSLois Curfman McInnes if (snes->type == SNES_NONLINEAR_EQUATIONS) { 29177c4ece6SBarry Smith PetscPrintf(snes->comm, 292d31a3109SLois Curfman McInnes " Options for solving systems of nonlinear equations only:\n"); 29377c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 29477c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 295ef1dfb11SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 29693c39befSBarry Smith PetscPrintf(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p); 29777c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 29877c4ece6SBarry Smith PetscPrintf(snes->comm, 299b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 30077c4ece6SBarry Smith PetscPrintf(snes->comm, 301b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 30277c4ece6SBarry Smith PetscPrintf(snes->comm, 303b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 30477c4ece6SBarry Smith PetscPrintf(snes->comm, 305b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 30677c4ece6SBarry Smith PetscPrintf(snes->comm, 307b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 30877c4ece6SBarry Smith PetscPrintf(snes->comm, 309b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 31077c4ece6SBarry Smith PetscPrintf(snes->comm, 311b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 312b18e04deSLois Curfman McInnes } 313b18e04deSLois Curfman McInnes else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 31477c4ece6SBarry Smith PetscPrintf(snes->comm, 315d31a3109SLois Curfman McInnes " Options for solving unconstrained minimization problems only:\n"); 316b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 31777c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 31877c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 319b18e04deSLois Curfman McInnes } 3204537a946SLois Curfman McInnes PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 32177c4ece6SBarry Smith PetscPrintf(snes->comm,"a particular method\n"); 3226daaf66cSBarry Smith if (snes->printhelp) (*snes->printhelp)(snes,p); 3239b94acceSBarry Smith return 0; 3249b94acceSBarry Smith } 325a847f771SSatish Balay 3265615d1e5SSatish Balay #undef __FUNC__ 327d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext" 3289b94acceSBarry Smith /*@ 3299b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 3309b94acceSBarry Smith the nonlinear solvers. 3319b94acceSBarry Smith 3329b94acceSBarry Smith Input Parameters: 3339b94acceSBarry Smith . snes - the SNES context 3349b94acceSBarry Smith . usrP - optional user context 3359b94acceSBarry Smith 3369b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3379b94acceSBarry Smith 3389b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3399b94acceSBarry Smith @*/ 3409b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3419b94acceSBarry Smith { 34277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3439b94acceSBarry Smith snes->user = usrP; 3449b94acceSBarry Smith return 0; 3459b94acceSBarry Smith } 34674679c65SBarry Smith 3475615d1e5SSatish Balay #undef __FUNC__ 348d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext" 3499b94acceSBarry Smith /*@C 3509b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3519b94acceSBarry Smith nonlinear solvers. 3529b94acceSBarry Smith 3539b94acceSBarry Smith Input Parameter: 3549b94acceSBarry Smith . snes - SNES context 3559b94acceSBarry Smith 3569b94acceSBarry Smith Output Parameter: 3579b94acceSBarry Smith . usrP - user context 3589b94acceSBarry Smith 3599b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3609b94acceSBarry Smith 3619b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3629b94acceSBarry Smith @*/ 3639b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3649b94acceSBarry Smith { 36577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3669b94acceSBarry Smith *usrP = snes->user; 3679b94acceSBarry Smith return 0; 3689b94acceSBarry Smith } 36974679c65SBarry Smith 3705615d1e5SSatish Balay #undef __FUNC__ 371d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber" 3729b94acceSBarry Smith /*@ 3739b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3749b94acceSBarry Smith nonlinear solver. 3759b94acceSBarry Smith 3769b94acceSBarry Smith Input Parameter: 3779b94acceSBarry Smith . snes - SNES context 3789b94acceSBarry Smith 3799b94acceSBarry Smith Output Parameter: 3809b94acceSBarry Smith . iter - iteration number 3819b94acceSBarry Smith 3829b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3839b94acceSBarry Smith @*/ 3849b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3859b94acceSBarry Smith { 38677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 38774679c65SBarry Smith PetscValidIntPointer(iter); 3889b94acceSBarry Smith *iter = snes->iter; 3899b94acceSBarry Smith return 0; 3909b94acceSBarry Smith } 39174679c65SBarry Smith 3925615d1e5SSatish Balay #undef __FUNC__ 3935615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 3949b94acceSBarry Smith /*@ 3959b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3969b94acceSBarry Smith with SNESSSetFunction(). 3979b94acceSBarry Smith 3989b94acceSBarry Smith Input Parameter: 3999b94acceSBarry Smith . snes - SNES context 4009b94acceSBarry Smith 4019b94acceSBarry Smith Output Parameter: 4029b94acceSBarry Smith . fnorm - 2-norm of function 4039b94acceSBarry Smith 4049b94acceSBarry Smith Note: 4059b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 406a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 407a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 4089b94acceSBarry Smith 4099b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 410a86d99e1SLois Curfman McInnes 411a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 4129b94acceSBarry Smith @*/ 4139b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 4149b94acceSBarry Smith { 41577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 41674679c65SBarry Smith PetscValidScalarPointer(fnorm); 41774679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 418d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 41974679c65SBarry Smith } 4209b94acceSBarry Smith *fnorm = snes->norm; 4219b94acceSBarry Smith return 0; 4229b94acceSBarry Smith } 42374679c65SBarry Smith 4245615d1e5SSatish Balay #undef __FUNC__ 4255615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 4269b94acceSBarry Smith /*@ 4279b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 4289b94acceSBarry Smith with SNESSSetGradient(). 4299b94acceSBarry Smith 4309b94acceSBarry Smith Input Parameter: 4319b94acceSBarry Smith . snes - SNES context 4329b94acceSBarry Smith 4339b94acceSBarry Smith Output Parameter: 4349b94acceSBarry Smith . fnorm - 2-norm of gradient 4359b94acceSBarry Smith 4369b94acceSBarry Smith Note: 4379b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 438a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 439a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4409b94acceSBarry Smith 4419b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 442a86d99e1SLois Curfman McInnes 443a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4449b94acceSBarry Smith @*/ 4459b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4469b94acceSBarry Smith { 44777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 44874679c65SBarry Smith PetscValidScalarPointer(gnorm); 44974679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 450d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 45174679c65SBarry Smith } 4529b94acceSBarry Smith *gnorm = snes->norm; 4539b94acceSBarry Smith return 0; 4549b94acceSBarry Smith } 45574679c65SBarry Smith 4565615d1e5SSatish Balay #undef __FUNC__ 457d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 4589b94acceSBarry Smith /*@ 4599b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4609b94acceSBarry Smith attempted by the nonlinear solver. 4619b94acceSBarry Smith 4629b94acceSBarry Smith Input Parameter: 4639b94acceSBarry Smith . snes - SNES context 4649b94acceSBarry Smith 4659b94acceSBarry Smith Output Parameter: 4669b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4679b94acceSBarry Smith 468c96a6f78SLois Curfman McInnes Notes: 469c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 470c96a6f78SLois Curfman McInnes 4719b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4729b94acceSBarry Smith @*/ 4739b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4749b94acceSBarry Smith { 47577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 47674679c65SBarry Smith PetscValidIntPointer(nfails); 4779b94acceSBarry Smith *nfails = snes->nfailures; 4789b94acceSBarry Smith return 0; 4799b94acceSBarry Smith } 480a847f771SSatish Balay 4815615d1e5SSatish Balay #undef __FUNC__ 482d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" 483c96a6f78SLois Curfman McInnes /*@ 484c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 485c96a6f78SLois Curfman McInnes used by the nonlinear solver. 486c96a6f78SLois Curfman McInnes 487c96a6f78SLois Curfman McInnes Input Parameter: 488c96a6f78SLois Curfman McInnes . snes - SNES context 489c96a6f78SLois Curfman McInnes 490c96a6f78SLois Curfman McInnes Output Parameter: 491c96a6f78SLois Curfman McInnes . lits - number of linear iterations 492c96a6f78SLois Curfman McInnes 493c96a6f78SLois Curfman McInnes Notes: 494c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 495c96a6f78SLois Curfman McInnes 496c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 497c96a6f78SLois Curfman McInnes @*/ 49886bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 499c96a6f78SLois Curfman McInnes { 500c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 501c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 502c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 503c96a6f78SLois Curfman McInnes return 0; 504c96a6f78SLois Curfman McInnes } 505c96a6f78SLois Curfman McInnes 506c96a6f78SLois Curfman McInnes #undef __FUNC__ 507d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES" 5089b94acceSBarry Smith /*@C 5099b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 5109b94acceSBarry Smith 5119b94acceSBarry Smith Input Parameter: 5129b94acceSBarry Smith . snes - the SNES context 5139b94acceSBarry Smith 5149b94acceSBarry Smith Output Parameter: 5159b94acceSBarry Smith . sles - the SLES context 5169b94acceSBarry Smith 5179b94acceSBarry Smith Notes: 5189b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 5199b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5209b94acceSBarry Smith KSP and PC contexts as well. 5219b94acceSBarry Smith 5229b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 5239b94acceSBarry Smith 5249b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 5259b94acceSBarry Smith @*/ 5269b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 5279b94acceSBarry Smith { 52877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5299b94acceSBarry Smith *sles = snes->sles; 5309b94acceSBarry Smith return 0; 5319b94acceSBarry Smith } 5329b94acceSBarry Smith /* -----------------------------------------------------------*/ 5335615d1e5SSatish Balay #undef __FUNC__ 5345615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 5359b94acceSBarry Smith /*@C 5369b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5379b94acceSBarry Smith 5389b94acceSBarry Smith Input Parameter: 5399b94acceSBarry Smith . comm - MPI communicator 5409b94acceSBarry Smith . type - type of method, one of 5419b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS 5429b94acceSBarry Smith $ (for systems of nonlinear equations) 5439b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION 5449b94acceSBarry Smith $ (for unconstrained minimization) 5459b94acceSBarry Smith 5469b94acceSBarry Smith Output Parameter: 5479b94acceSBarry Smith . outsnes - the new SNES context 5489b94acceSBarry Smith 549c1f60f51SBarry Smith Options Database Key: 55019bd6747SLois Curfman McInnes $ -snes_mf - use default matrix-free Jacobian-vector products, 55119bd6747SLois Curfman McInnes $ and no preconditioning matrix 55219bd6747SLois Curfman McInnes $ -snes_mf_operator - use default matrix-free Jacobian-vector 55319bd6747SLois Curfman McInnes $ products, and a user-provided preconditioning matrix 55419bd6747SLois Curfman McInnes $ as set by SNESSetJacobian() 55519bd6747SLois Curfman McInnes $ -snes_fd - use (slow!) finite differences to compute Jacobian 556c1f60f51SBarry Smith 5579b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5589b94acceSBarry Smith 55963a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 5609b94acceSBarry Smith @*/ 5614b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 5629b94acceSBarry Smith { 5639b94acceSBarry Smith int ierr; 5649b94acceSBarry Smith SNES snes; 5659b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 56637fcc0dbSBarry Smith 5679b94acceSBarry Smith *outsnes = 0; 5682a0acf01SLois Curfman McInnes if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 569d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 570d4bb536fSBarry Smith PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm,SNESDestroy,SNESView); 5719b94acceSBarry Smith PLogObjectCreate(snes); 5729b94acceSBarry Smith snes->max_its = 50; 5739750a799SBarry Smith snes->max_funcs = 10000; 5749b94acceSBarry Smith snes->norm = 0.0; 5755a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 576b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 577b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5789b94acceSBarry Smith snes->atol = 1.e-10; 5795a2d0531SBarry Smith } 580b4874afaSBarry Smith else { 581b4874afaSBarry Smith snes->rtol = 1.e-8; 582b4874afaSBarry Smith snes->ttol = 0.0; 583b4874afaSBarry Smith snes->atol = 1.e-50; 584b4874afaSBarry Smith } 5859b94acceSBarry Smith snes->xtol = 1.e-8; 586b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5879b94acceSBarry Smith snes->nfuncs = 0; 5889b94acceSBarry Smith snes->nfailures = 0; 5897a00f4a9SLois Curfman McInnes snes->linear_its = 0; 590639f9d9dSBarry Smith snes->numbermonitors = 0; 5919b94acceSBarry Smith snes->data = 0; 5929b94acceSBarry Smith snes->view = 0; 5939b94acceSBarry Smith snes->computeumfunction = 0; 5949b94acceSBarry Smith snes->umfunP = 0; 5959b94acceSBarry Smith snes->fc = 0; 5969b94acceSBarry Smith snes->deltatol = 1.e-12; 5979b94acceSBarry Smith snes->fmin = -1.e30; 5989b94acceSBarry Smith snes->method_class = type; 5999b94acceSBarry Smith snes->set_method_called = 0; 600a703fe33SLois Curfman McInnes snes->setup_called = 0; 6019b94acceSBarry Smith snes->ksp_ewconv = 0; 6027d1a2b2bSBarry Smith snes->type = -1; 6036f24a144SLois Curfman McInnes snes->xmonitor = 0; 6046f24a144SLois Curfman McInnes snes->vwork = 0; 6056f24a144SLois Curfman McInnes snes->nwork = 0; 606c9005455SLois Curfman McInnes snes->conv_hist_size = 0; 607c9005455SLois Curfman McInnes snes->conv_act_size = 0; 608c9005455SLois Curfman McInnes snes->conv_hist = 0; 6099b94acceSBarry Smith 6109b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 6110452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 612eed86810SBarry Smith PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 6139b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6149b94acceSBarry Smith kctx->version = 2; 6159b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6169b94acceSBarry Smith this was too large for some test cases */ 6179b94acceSBarry Smith kctx->rtol_last = 0; 6189b94acceSBarry Smith kctx->rtol_max = .9; 6199b94acceSBarry Smith kctx->gamma = 1.0; 6209b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6219b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6229b94acceSBarry Smith kctx->threshold = .1; 6239b94acceSBarry Smith kctx->lresid_last = 0; 6249b94acceSBarry Smith kctx->norm_last = 0; 6259b94acceSBarry Smith 6269b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 6279b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 6285334005bSBarry Smith 6299b94acceSBarry Smith *outsnes = snes; 6309b94acceSBarry Smith return 0; 6319b94acceSBarry Smith } 6329b94acceSBarry Smith 6339b94acceSBarry Smith /* --------------------------------------------------------------- */ 6345615d1e5SSatish Balay #undef __FUNC__ 635d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction" 6369b94acceSBarry Smith /*@C 6379b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6389b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6399b94acceSBarry Smith equations. 6409b94acceSBarry Smith 6419b94acceSBarry Smith Input Parameters: 6429b94acceSBarry Smith . snes - the SNES context 6439b94acceSBarry Smith . func - function evaluation routine 6449b94acceSBarry Smith . r - vector to store function value 6452cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 6462cd2dfdcSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6479b94acceSBarry Smith 6489b94acceSBarry Smith Calling sequence of func: 649313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Vec f,void *ctx); 6509b94acceSBarry Smith 6519b94acceSBarry Smith . x - input vector 652313e4042SLois Curfman McInnes . f - function vector 6532cd2dfdcSLois Curfman McInnes . ctx - optional user-defined function context 6549b94acceSBarry Smith 6559b94acceSBarry Smith Notes: 6569b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6579b94acceSBarry Smith $ f'(x) x = -f(x), 6589b94acceSBarry Smith $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 6599b94acceSBarry Smith 6609b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6619b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6629b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 6639b94acceSBarry Smith 6649b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6659b94acceSBarry Smith 666a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6679b94acceSBarry Smith @*/ 6685334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 6699b94acceSBarry Smith { 67077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 671e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6729b94acceSBarry Smith snes->computefunction = func; 6739b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6749b94acceSBarry Smith snes->funP = ctx; 6759b94acceSBarry Smith return 0; 6769b94acceSBarry Smith } 6779b94acceSBarry Smith 6785615d1e5SSatish Balay #undef __FUNC__ 6795615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 6809b94acceSBarry Smith /*@ 6819b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6829b94acceSBarry Smith SNESSetFunction(). 6839b94acceSBarry Smith 6849b94acceSBarry Smith Input Parameters: 6859b94acceSBarry Smith . snes - the SNES context 6869b94acceSBarry Smith . x - input vector 6879b94acceSBarry Smith 6889b94acceSBarry Smith Output Parameter: 6893638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6909b94acceSBarry Smith 6911bffabb2SLois Curfman McInnes Notes: 6929b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6939b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6949b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6959b94acceSBarry Smith 6969b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6979b94acceSBarry Smith 698a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6999b94acceSBarry Smith @*/ 7009b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 7019b94acceSBarry Smith { 7029b94acceSBarry Smith int ierr; 7039b94acceSBarry Smith 704d4bb536fSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 705e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 706d4bb536fSBarry Smith } 7079b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 7089b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 709ae3c334cSLois Curfman McInnes snes->nfuncs++; 7109b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 7119b94acceSBarry Smith return 0; 7129b94acceSBarry Smith } 7139b94acceSBarry Smith 7145615d1e5SSatish Balay #undef __FUNC__ 715d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" 7169b94acceSBarry Smith /*@C 7179b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 7189b94acceSBarry Smith unconstrained minimization. 7199b94acceSBarry Smith 7209b94acceSBarry Smith Input Parameters: 7219b94acceSBarry Smith . snes - the SNES context 7229b94acceSBarry Smith . func - function evaluation routine 723044dda88SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 724044dda88SLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 7259b94acceSBarry Smith 7269b94acceSBarry Smith Calling sequence of func: 7279b94acceSBarry Smith . func (SNES snes,Vec x,double *f,void *ctx); 7289b94acceSBarry Smith 7299b94acceSBarry Smith . x - input vector 7309b94acceSBarry Smith . f - function 731044dda88SLois Curfman McInnes . ctx - [optional] user-defined function context 7329b94acceSBarry Smith 7339b94acceSBarry Smith Notes: 7349b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7359b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7369b94acceSBarry Smith SNESSetFunction(). 7379b94acceSBarry Smith 7389b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 7399b94acceSBarry Smith 740a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 741a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 7429b94acceSBarry Smith @*/ 7439b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 7449b94acceSBarry Smith void *ctx) 7459b94acceSBarry Smith { 74677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 747e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7489b94acceSBarry Smith snes->computeumfunction = func; 7499b94acceSBarry Smith snes->umfunP = ctx; 7509b94acceSBarry Smith return 0; 7519b94acceSBarry Smith } 7529b94acceSBarry Smith 7535615d1e5SSatish Balay #undef __FUNC__ 7545615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 7559b94acceSBarry Smith /*@ 7569b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 7579b94acceSBarry Smith set with SNESSetMinimizationFunction(). 7589b94acceSBarry Smith 7599b94acceSBarry Smith Input Parameters: 7609b94acceSBarry Smith . snes - the SNES context 7619b94acceSBarry Smith . x - input vector 7629b94acceSBarry Smith 7639b94acceSBarry Smith Output Parameter: 7649b94acceSBarry Smith . y - function value 7659b94acceSBarry Smith 7669b94acceSBarry Smith Notes: 7679b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 7689b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7699b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 770a86d99e1SLois Curfman McInnes 771a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 772a86d99e1SLois Curfman McInnes 773a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 774a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 7759b94acceSBarry Smith @*/ 7769b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 7779b94acceSBarry Smith { 7789b94acceSBarry Smith int ierr; 779e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7809b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 7819b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 782ae3c334cSLois Curfman McInnes snes->nfuncs++; 7839b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7849b94acceSBarry Smith return 0; 7859b94acceSBarry Smith } 7869b94acceSBarry Smith 7875615d1e5SSatish Balay #undef __FUNC__ 788d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient" 7899b94acceSBarry Smith /*@C 7909b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 7919b94acceSBarry Smith vector for use by the SNES routines. 7929b94acceSBarry Smith 7939b94acceSBarry Smith Input Parameters: 7949b94acceSBarry Smith . snes - the SNES context 7959b94acceSBarry Smith . func - function evaluation routine 796044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 797044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 7989b94acceSBarry Smith . r - vector to store gradient value 7999b94acceSBarry Smith 8009b94acceSBarry Smith Calling sequence of func: 8019b94acceSBarry Smith . func (SNES, Vec x, Vec g, void *ctx); 8029b94acceSBarry Smith 8039b94acceSBarry Smith . x - input vector 8049b94acceSBarry Smith . g - gradient vector 805044dda88SLois Curfman McInnes . ctx - optional user-defined gradient context 8069b94acceSBarry Smith 8079b94acceSBarry Smith Notes: 8089b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 8099b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 8109b94acceSBarry Smith SNESSetFunction(). 8119b94acceSBarry Smith 8129b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 8139b94acceSBarry Smith 814a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 815a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 8169b94acceSBarry Smith @*/ 81774679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 8189b94acceSBarry Smith { 81977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 820e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8219b94acceSBarry Smith snes->computefunction = func; 8229b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 8239b94acceSBarry Smith snes->funP = ctx; 8249b94acceSBarry Smith return 0; 8259b94acceSBarry Smith } 8269b94acceSBarry Smith 8275615d1e5SSatish Balay #undef __FUNC__ 8285615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 8299b94acceSBarry Smith /*@ 830a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 831a86d99e1SLois Curfman McInnes SNESSetGradient(). 8329b94acceSBarry Smith 8339b94acceSBarry Smith Input Parameters: 8349b94acceSBarry Smith . snes - the SNES context 8359b94acceSBarry Smith . x - input vector 8369b94acceSBarry Smith 8379b94acceSBarry Smith Output Parameter: 8389b94acceSBarry Smith . y - gradient vector 8399b94acceSBarry Smith 8409b94acceSBarry Smith Notes: 8419b94acceSBarry Smith SNESComputeGradient() is valid only for 8429b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 8439b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 844a86d99e1SLois Curfman McInnes 845a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 846a86d99e1SLois Curfman McInnes 847a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 848a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 8499b94acceSBarry Smith @*/ 8509b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 8519b94acceSBarry Smith { 8529b94acceSBarry Smith int ierr; 85374679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 854e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8559b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 8569b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 8579b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 8589b94acceSBarry Smith return 0; 8599b94acceSBarry Smith } 8609b94acceSBarry Smith 8615615d1e5SSatish Balay #undef __FUNC__ 8625615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 86362fef451SLois Curfman McInnes /*@ 86462fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 86562fef451SLois Curfman McInnes set with SNESSetJacobian(). 86662fef451SLois Curfman McInnes 86762fef451SLois Curfman McInnes Input Parameters: 86862fef451SLois Curfman McInnes . snes - the SNES context 86962fef451SLois Curfman McInnes . x - input vector 87062fef451SLois Curfman McInnes 87162fef451SLois Curfman McInnes Output Parameters: 87262fef451SLois Curfman McInnes . A - Jacobian matrix 87362fef451SLois Curfman McInnes . B - optional preconditioning matrix 87462fef451SLois Curfman McInnes . flag - flag indicating matrix structure 87562fef451SLois Curfman McInnes 87662fef451SLois Curfman McInnes Notes: 87762fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 87862fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 87962fef451SLois Curfman McInnes 880dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 881dc5a77f8SLois Curfman McInnes flag parameter. 88262fef451SLois Curfman McInnes 88362fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 88462fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 885005c665bSBarry Smith methods is SNESComputeHessian(). 88662fef451SLois Curfman McInnes 88762fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 88862fef451SLois Curfman McInnes 88962fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 89062fef451SLois Curfman McInnes @*/ 8919b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8929b94acceSBarry Smith { 8939b94acceSBarry Smith int ierr; 89474679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 895e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 8969b94acceSBarry Smith if (!snes->computejacobian) return 0; 8979b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 898c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 8999b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 9009b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 9016d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 90277c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 90377c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9049b94acceSBarry Smith return 0; 9059b94acceSBarry Smith } 9069b94acceSBarry Smith 9075615d1e5SSatish Balay #undef __FUNC__ 9085615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 90962fef451SLois Curfman McInnes /*@ 91062fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 91162fef451SLois Curfman McInnes set with SNESSetHessian(). 91262fef451SLois Curfman McInnes 91362fef451SLois Curfman McInnes Input Parameters: 91462fef451SLois Curfman McInnes . snes - the SNES context 91562fef451SLois Curfman McInnes . x - input vector 91662fef451SLois Curfman McInnes 91762fef451SLois Curfman McInnes Output Parameters: 91862fef451SLois Curfman McInnes . A - Hessian matrix 91962fef451SLois Curfman McInnes . B - optional preconditioning matrix 92062fef451SLois Curfman McInnes . flag - flag indicating matrix structure 92162fef451SLois Curfman McInnes 92262fef451SLois Curfman McInnes Notes: 92362fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 92462fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 92562fef451SLois Curfman McInnes 926dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 927dc5a77f8SLois Curfman McInnes flag parameter. 92862fef451SLois Curfman McInnes 92962fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 93062fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 93162fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 93262fef451SLois Curfman McInnes 93362fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 93462fef451SLois Curfman McInnes 935a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 936a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 93762fef451SLois Curfman McInnes @*/ 93862fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 9399b94acceSBarry Smith { 9409b94acceSBarry Smith int ierr; 94174679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 942e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 9439b94acceSBarry Smith if (!snes->computejacobian) return 0; 94462fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 9450f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 94662fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 94762fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 9480f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 94977c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 95077c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9519b94acceSBarry Smith return 0; 9529b94acceSBarry Smith } 9539b94acceSBarry Smith 9545615d1e5SSatish Balay #undef __FUNC__ 955d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian" 9569b94acceSBarry Smith /*@C 9579b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 958044dda88SLois Curfman McInnes location to store the matrix. 9599b94acceSBarry Smith 9609b94acceSBarry Smith Input Parameters: 9619b94acceSBarry Smith . snes - the SNES context 9629b94acceSBarry Smith . A - Jacobian matrix 9639b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 9649b94acceSBarry Smith . func - Jacobian evaluation routine 9652cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 9662cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 9679b94acceSBarry Smith 9689b94acceSBarry Smith Calling sequence of func: 969313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9709b94acceSBarry Smith 9719b94acceSBarry Smith . x - input vector 9729b94acceSBarry Smith . A - Jacobian matrix 9739b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 974ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 975ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 9762cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Jacobian context 9779b94acceSBarry Smith 9789b94acceSBarry Smith Notes: 979dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 9802cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 981ac21db08SLois Curfman McInnes 982ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9839b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9849b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9859b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9869b94acceSBarry Smith throughout the global iterations. 9879b94acceSBarry Smith 9889b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9899b94acceSBarry Smith 990ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 9919b94acceSBarry Smith @*/ 9929b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 9939b94acceSBarry Smith MatStructure*,void*),void *ctx) 9949b94acceSBarry Smith { 99577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 996d4bb536fSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 9979b94acceSBarry Smith snes->computejacobian = func; 9989b94acceSBarry Smith snes->jacP = ctx; 9999b94acceSBarry Smith snes->jacobian = A; 10009b94acceSBarry Smith snes->jacobian_pre = B; 10019b94acceSBarry Smith return 0; 10029b94acceSBarry Smith } 100362fef451SLois Curfman McInnes 10045615d1e5SSatish Balay #undef __FUNC__ 1005d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian" 1006b4fd4287SBarry Smith /*@ 1007b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1008b4fd4287SBarry Smith provided context for evaluating the Jacobian. 1009b4fd4287SBarry Smith 1010b4fd4287SBarry Smith Input Parameter: 1011b4fd4287SBarry Smith . snes - the nonlinear solver context 1012b4fd4287SBarry Smith 1013b4fd4287SBarry Smith Output Parameters: 1014b4fd4287SBarry Smith . A - location to stash Jacobian matrix (or PETSC_NULL) 1015b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 1016b4fd4287SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 1017b4fd4287SBarry Smith 1018b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 1019b4fd4287SBarry Smith @*/ 1020b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1021b4fd4287SBarry Smith { 102274679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 1023e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 1024b4fd4287SBarry Smith if (A) *A = snes->jacobian; 1025b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 1026b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 1027b4fd4287SBarry Smith return 0; 1028b4fd4287SBarry Smith } 1029b4fd4287SBarry Smith 10305615d1e5SSatish Balay #undef __FUNC__ 1031d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian" 10329b94acceSBarry Smith /*@C 10339b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 1034044dda88SLois Curfman McInnes location to store the matrix. 10359b94acceSBarry Smith 10369b94acceSBarry Smith Input Parameters: 10379b94acceSBarry Smith . snes - the SNES context 10389b94acceSBarry Smith . A - Hessian matrix 10399b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 10409b94acceSBarry Smith . func - Jacobian evaluation routine 1041313e4042SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 1042313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 10439b94acceSBarry Smith 10449b94acceSBarry Smith Calling sequence of func: 1045313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 10469b94acceSBarry Smith 10479b94acceSBarry Smith . x - input vector 10489b94acceSBarry Smith . A - Hessian matrix 10499b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1050ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1051ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 10522cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Hessian context 10539b94acceSBarry Smith 10549b94acceSBarry Smith Notes: 1055dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10562cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1057ac21db08SLois Curfman McInnes 10589b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 10599b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 10609b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10619b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10629b94acceSBarry Smith throughout the global iterations. 10639b94acceSBarry Smith 10649b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 10659b94acceSBarry Smith 1066ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 10679b94acceSBarry Smith @*/ 10689b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 10699b94acceSBarry Smith MatStructure*,void*),void *ctx) 10709b94acceSBarry Smith { 107177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1072d4bb536fSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1073e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1074d4bb536fSBarry Smith } 10759b94acceSBarry Smith snes->computejacobian = func; 10769b94acceSBarry Smith snes->jacP = ctx; 10779b94acceSBarry Smith snes->jacobian = A; 10789b94acceSBarry Smith snes->jacobian_pre = B; 10799b94acceSBarry Smith return 0; 10809b94acceSBarry Smith } 10819b94acceSBarry Smith 10825615d1e5SSatish Balay #undef __FUNC__ 1083d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian" 108462fef451SLois Curfman McInnes /*@ 108562fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 108662fef451SLois Curfman McInnes provided context for evaluating the Hessian. 108762fef451SLois Curfman McInnes 108862fef451SLois Curfman McInnes Input Parameter: 108962fef451SLois Curfman McInnes . snes - the nonlinear solver context 109062fef451SLois Curfman McInnes 109162fef451SLois Curfman McInnes Output Parameters: 109262fef451SLois Curfman McInnes . A - location to stash Hessian matrix (or PETSC_NULL) 109362fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 109462fef451SLois Curfman McInnes . ctx - location to stash Hessian ctx (or PETSC_NULL) 109562fef451SLois Curfman McInnes 109662fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 109762fef451SLois Curfman McInnes @*/ 109862fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 109962fef451SLois Curfman McInnes { 110074679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1101e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 110262fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 110362fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 110462fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 110562fef451SLois Curfman McInnes return 0; 110662fef451SLois Curfman McInnes } 110762fef451SLois Curfman McInnes 11089b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 11099b94acceSBarry Smith 11105615d1e5SSatish Balay #undef __FUNC__ 11115615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 11129b94acceSBarry Smith /*@ 11139b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1114272ac6f2SLois Curfman McInnes of a nonlinear solver. 11159b94acceSBarry Smith 11169b94acceSBarry Smith Input Parameter: 11179b94acceSBarry Smith . snes - the SNES context 11188ddd3da0SLois Curfman McInnes . x - the solution vector 11199b94acceSBarry Smith 1120272ac6f2SLois Curfman McInnes Notes: 1121272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1122272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1123272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1124272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1125272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1126272ac6f2SLois Curfman McInnes 11279b94acceSBarry Smith .keywords: SNES, nonlinear, setup 11289b94acceSBarry Smith 11299b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 11309b94acceSBarry Smith @*/ 11318ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 11329b94acceSBarry Smith { 1133272ac6f2SLois Curfman McInnes int ierr, flg; 113477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 113577c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 11368ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 11379b94acceSBarry Smith 1138c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1139c1f60f51SBarry Smith /* 1140c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1141dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1142c1f60f51SBarry Smith */ 1143c1f60f51SBarry Smith if (flg) { 1144c1f60f51SBarry Smith Mat J; 1145c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1146c1f60f51SBarry Smith PLogObjectParent(snes,J); 1147c1f60f51SBarry Smith snes->mfshell = J; 1148c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1149c1f60f51SBarry Smith snes->jacobian = J; 1150c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1151c1f60f51SBarry Smith } 1152c1f60f51SBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1153c1f60f51SBarry Smith snes->jacobian = J; 1154c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1155d4bb536fSBarry Smith } else { 1156d4bb536fSBarry Smith SETERRQ(1,0,"Method class doesn't support matrix-free operator option"); 1157d4bb536fSBarry Smith } 1158c1f60f51SBarry Smith } 1159272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1160c1f60f51SBarry Smith /* 1161dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1162c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1163c1f60f51SBarry Smith */ 116431615d04SLois Curfman McInnes if (flg) { 1165272ac6f2SLois Curfman McInnes Mat J; 1166272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1167272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1168272ac6f2SLois Curfman McInnes snes->mfshell = J; 116983e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 117083e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 117194a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 117283e56afdSLois Curfman McInnes } 117383e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 117483e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 117594a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1176d4bb536fSBarry Smith } else { 1177d4bb536fSBarry Smith SETERRQ(1,0,"Method class doesn't support matrix-free option"); 1178d4bb536fSBarry Smith } 1179272ac6f2SLois Curfman McInnes } 11809b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1181e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1182e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1183e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first"); 1184e3372554SBarry Smith if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector"); 1185a703fe33SLois Curfman McInnes 1186a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 118740191667SLois Curfman McInnes if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1188a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1189a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1190a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1191a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1192a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1193a703fe33SLois Curfman McInnes } 11949b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1195e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first"); 1196e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first"); 1197d4bb536fSBarry Smith if (!snes->computeumfunction) SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first"); 1198e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first"); 1199d4bb536fSBarry Smith } else { 1200d4bb536fSBarry Smith SETERRQ(1,0,"Unknown method class"); 1201d4bb536fSBarry Smith } 1202a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1203a703fe33SLois Curfman McInnes snes->setup_called = 1; 1204a703fe33SLois Curfman McInnes return 0; 12059b94acceSBarry Smith } 12069b94acceSBarry Smith 12075615d1e5SSatish Balay #undef __FUNC__ 1208d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy" 12099b94acceSBarry Smith /*@C 12109b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 12119b94acceSBarry Smith with SNESCreate(). 12129b94acceSBarry Smith 12139b94acceSBarry Smith Input Parameter: 12149b94acceSBarry Smith . snes - the SNES context 12159b94acceSBarry Smith 12169b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 12179b94acceSBarry Smith 121863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 12199b94acceSBarry Smith @*/ 12209b94acceSBarry Smith int SNESDestroy(SNES snes) 12219b94acceSBarry Smith { 12229b94acceSBarry Smith int ierr; 122377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1224d4bb536fSBarry Smith if (--snes->refct > 0) return 0; 1225d4bb536fSBarry Smith 12269750a799SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);} 12270452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 12289b94acceSBarry Smith if (snes->mfshell) MatDestroy(snes->mfshell); 12299b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 12303e2e452bSBarry Smith if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 12316f24a144SLois Curfman McInnes if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 12329b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 12330452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 12349b94acceSBarry Smith return 0; 12359b94acceSBarry Smith } 12369b94acceSBarry Smith 12379b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 12389b94acceSBarry Smith 12395615d1e5SSatish Balay #undef __FUNC__ 12405615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 12419b94acceSBarry Smith /*@ 1242d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 12439b94acceSBarry Smith 12449b94acceSBarry Smith Input Parameters: 12459b94acceSBarry Smith . snes - the SNES context 124633174efeSLois Curfman McInnes . atol - absolute convergence tolerance 124733174efeSLois Curfman McInnes . rtol - relative convergence tolerance 124833174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 124933174efeSLois Curfman McInnes of the change in the solution between steps 125033174efeSLois Curfman McInnes . maxit - maximum number of iterations 125133174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 12529b94acceSBarry Smith 125333174efeSLois Curfman McInnes Options Database Keys: 125433174efeSLois Curfman McInnes $ -snes_atol <atol> 125533174efeSLois Curfman McInnes $ -snes_rtol <rtol> 125633174efeSLois Curfman McInnes $ -snes_stol <stol> 125733174efeSLois Curfman McInnes $ -snes_max_it <maxit> 125833174efeSLois Curfman McInnes $ -snes_max_funcs <maxf> 12599b94acceSBarry Smith 1260d7a720efSLois Curfman McInnes Notes: 12619b94acceSBarry Smith The default maximum number of iterations is 50. 12629b94acceSBarry Smith The default maximum number of function evaluations is 1000. 12639b94acceSBarry Smith 126433174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 12659b94acceSBarry Smith 1266d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 12679b94acceSBarry Smith @*/ 1268d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 12699b94acceSBarry Smith { 127077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1271d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1272d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1273d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1274d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1275d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 12769b94acceSBarry Smith return 0; 12779b94acceSBarry Smith } 12789b94acceSBarry Smith 12795615d1e5SSatish Balay #undef __FUNC__ 12805615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 12819b94acceSBarry Smith /*@ 128233174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 128333174efeSLois Curfman McInnes 128433174efeSLois Curfman McInnes Input Parameters: 128533174efeSLois Curfman McInnes . snes - the SNES context 128633174efeSLois Curfman McInnes . atol - absolute convergence tolerance 128733174efeSLois Curfman McInnes . rtol - relative convergence tolerance 128833174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 128933174efeSLois Curfman McInnes of the change in the solution between steps 129033174efeSLois Curfman McInnes . maxit - maximum number of iterations 129133174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 129233174efeSLois Curfman McInnes 129333174efeSLois Curfman McInnes Notes: 129433174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 129533174efeSLois Curfman McInnes 129633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 129733174efeSLois Curfman McInnes 129833174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 129933174efeSLois Curfman McInnes @*/ 130033174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 130133174efeSLois Curfman McInnes { 130233174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 130333174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 130433174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 130533174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 130633174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 130733174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 130833174efeSLois Curfman McInnes return 0; 130933174efeSLois Curfman McInnes } 131033174efeSLois Curfman McInnes 13115615d1e5SSatish Balay #undef __FUNC__ 13125615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 131333174efeSLois Curfman McInnes /*@ 13149b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 13159b94acceSBarry Smith 13169b94acceSBarry Smith Input Parameters: 13179b94acceSBarry Smith . snes - the SNES context 13189b94acceSBarry Smith . tol - tolerance 13199b94acceSBarry Smith 13209b94acceSBarry Smith Options Database Key: 1321d7a720efSLois Curfman McInnes $ -snes_trtol <tol> 13229b94acceSBarry Smith 13239b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 13249b94acceSBarry Smith 1325d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 13269b94acceSBarry Smith @*/ 13279b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 13289b94acceSBarry Smith { 132977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13309b94acceSBarry Smith snes->deltatol = tol; 13319b94acceSBarry Smith return 0; 13329b94acceSBarry Smith } 13339b94acceSBarry Smith 13345615d1e5SSatish Balay #undef __FUNC__ 13355615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 13369b94acceSBarry Smith /*@ 133777c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 13389b94acceSBarry Smith for unconstrained minimization solvers. 13399b94acceSBarry Smith 13409b94acceSBarry Smith Input Parameters: 13419b94acceSBarry Smith . snes - the SNES context 13429b94acceSBarry Smith . ftol - minimum function tolerance 13439b94acceSBarry Smith 13449b94acceSBarry Smith Options Database Key: 1345d7a720efSLois Curfman McInnes $ -snes_fmin <ftol> 13469b94acceSBarry Smith 13479b94acceSBarry Smith Note: 134877c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 13499b94acceSBarry Smith methods only. 13509b94acceSBarry Smith 13519b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 13529b94acceSBarry Smith 1353d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 13549b94acceSBarry Smith @*/ 135577c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 13569b94acceSBarry Smith { 135777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13589b94acceSBarry Smith snes->fmin = ftol; 13599b94acceSBarry Smith return 0; 13609b94acceSBarry Smith } 13619b94acceSBarry Smith 13629b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 13639b94acceSBarry Smith 13645615d1e5SSatish Balay #undef __FUNC__ 1365d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor" 13669b94acceSBarry Smith /*@C 13679b94acceSBarry Smith SNESSetMonitor - Sets the function that is to be used at every 13689b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 13699b94acceSBarry Smith progress. 13709b94acceSBarry Smith 13719b94acceSBarry Smith Input Parameters: 13729b94acceSBarry Smith . snes - the SNES context 13739b94acceSBarry Smith . func - monitoring routine 1374044dda88SLois Curfman McInnes . mctx - [optional] user-defined context for private data for the 1375044dda88SLois Curfman McInnes monitor routine (may be PETSC_NULL) 13769b94acceSBarry Smith 13779b94acceSBarry Smith Calling sequence of func: 1378682d7d0cSBarry Smith int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 13799b94acceSBarry Smith 13809b94acceSBarry Smith $ snes - the SNES context 13819b94acceSBarry Smith $ its - iteration number 1382044dda88SLois Curfman McInnes $ mctx - [optional] monitoring context 13839b94acceSBarry Smith $ 13849b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13859b94acceSBarry Smith $ norm - 2-norm function value (may be estimated) 13869b94acceSBarry Smith $ 13879b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13889b94acceSBarry Smith $ norm - 2-norm gradient value (may be estimated) 13899b94acceSBarry Smith 13909665c990SLois Curfman McInnes Options Database Keys: 13919665c990SLois Curfman McInnes $ -snes_monitor : sets SNESDefaultMonitor() 13929665c990SLois Curfman McInnes $ -snes_xmonitor : sets line graph monitor, 13939665c990SLois Curfman McInnes $ uses SNESLGMonitorCreate() 13949665c990SLois Curfman McInnes $ -snes_cancelmonitors : cancels all monitors that have 13959665c990SLois Curfman McInnes $ been hardwired into a code by 13969665c990SLois Curfman McInnes $ calls to SNESSetMonitor(), but 13979665c990SLois Curfman McInnes $ does not cancel those set via 13989665c990SLois Curfman McInnes $ the options database. 13999665c990SLois Curfman McInnes 14009665c990SLois Curfman McInnes 1401639f9d9dSBarry Smith Notes: 14026bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 14036bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 14046bc08f3fSLois Curfman McInnes order in which they were set. 1405639f9d9dSBarry Smith 14069b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 14079b94acceSBarry Smith 14089b94acceSBarry Smith .seealso: SNESDefaultMonitor() 14099b94acceSBarry Smith @*/ 141074679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 14119b94acceSBarry Smith { 1412639f9d9dSBarry Smith if (!func) { 1413639f9d9dSBarry Smith snes->numbermonitors = 0; 1414639f9d9dSBarry Smith return 0; 1415639f9d9dSBarry Smith } 1416639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1417e3372554SBarry Smith SETERRQ(1,0,"Too many monitors set"); 1418639f9d9dSBarry Smith } 1419639f9d9dSBarry Smith 1420639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1421639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 14229b94acceSBarry Smith return 0; 14239b94acceSBarry Smith } 14249b94acceSBarry Smith 14255615d1e5SSatish Balay #undef __FUNC__ 1426d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest" 14279b94acceSBarry Smith /*@C 14289b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 14299b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 14309b94acceSBarry Smith 14319b94acceSBarry Smith Input Parameters: 14329b94acceSBarry Smith . snes - the SNES context 14339b94acceSBarry Smith . func - routine to test for convergence 1434044dda88SLois Curfman McInnes . cctx - [optional] context for private data for the convergence routine 1435044dda88SLois Curfman McInnes (may be PETSC_NULL) 14369b94acceSBarry Smith 14379b94acceSBarry Smith Calling sequence of func: 14389b94acceSBarry Smith int func (SNES snes,double xnorm,double gnorm, 14399b94acceSBarry Smith double f,void *cctx) 14409b94acceSBarry Smith 14419b94acceSBarry Smith $ snes - the SNES context 1442044dda88SLois Curfman McInnes $ cctx - [optional] convergence context 14439b94acceSBarry Smith $ xnorm - 2-norm of current iterate 14449b94acceSBarry Smith $ 14459b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 14469b94acceSBarry Smith $ gnorm - 2-norm of current step 14479b94acceSBarry Smith $ f - 2-norm of function 14489b94acceSBarry Smith $ 14499b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 14509b94acceSBarry Smith $ gnorm - 2-norm of current gradient 14519b94acceSBarry Smith $ f - function value 14529b94acceSBarry Smith 14539b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 14549b94acceSBarry Smith 145540191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 145640191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 14579b94acceSBarry Smith @*/ 145874679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 14599b94acceSBarry Smith { 14609b94acceSBarry Smith (snes)->converged = func; 14619b94acceSBarry Smith (snes)->cnvP = cctx; 14629b94acceSBarry Smith return 0; 14639b94acceSBarry Smith } 14649b94acceSBarry Smith 14655615d1e5SSatish Balay #undef __FUNC__ 1466d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" 1467c9005455SLois Curfman McInnes /*@ 1468c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1469c9005455SLois Curfman McInnes 1470c9005455SLois Curfman McInnes Input Parameters: 1471c9005455SLois Curfman McInnes . snes - iterative context obtained from SNESCreate() 1472c9005455SLois Curfman McInnes . a - array to hold history 1473c9005455SLois Curfman McInnes . na - size of a 1474c9005455SLois Curfman McInnes 1475c9005455SLois Curfman McInnes Notes: 1476c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1477c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1478c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1479c9005455SLois Curfman McInnes at each step. 1480c9005455SLois Curfman McInnes 1481c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1482c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1483c9005455SLois Curfman McInnes during the section of code that is being timed. 1484c9005455SLois Curfman McInnes 1485c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1486c9005455SLois Curfman McInnes @*/ 1487c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na) 1488c9005455SLois Curfman McInnes { 1489c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1490c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1491c9005455SLois Curfman McInnes snes->conv_hist = a; 1492c9005455SLois Curfman McInnes snes->conv_hist_size = na; 1493c9005455SLois Curfman McInnes return 0; 1494c9005455SLois Curfman McInnes } 1495c9005455SLois Curfman McInnes 1496c9005455SLois Curfman McInnes #undef __FUNC__ 14975615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 14989b94acceSBarry Smith /* 14999b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 15009b94acceSBarry Smith positive parameter delta. 15019b94acceSBarry Smith 15029b94acceSBarry Smith Input Parameters: 15039b94acceSBarry Smith . snes - the SNES context 15049b94acceSBarry Smith . y - approximate solution of linear system 15059b94acceSBarry Smith . fnorm - 2-norm of current function 15069b94acceSBarry Smith . delta - trust region size 15079b94acceSBarry Smith 15089b94acceSBarry Smith Output Parameters: 15099b94acceSBarry Smith . gpnorm - predicted function norm at the new point, assuming local 15109b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 15119b94acceSBarry Smith region, and exceeds zero otherwise. 15129b94acceSBarry Smith . ynorm - 2-norm of the step 15139b94acceSBarry Smith 15149b94acceSBarry Smith Note: 151540191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 15169b94acceSBarry Smith is set to be the maximum allowable step size. 15179b94acceSBarry Smith 15189b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 15199b94acceSBarry Smith */ 15209b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 15219b94acceSBarry Smith double *gpnorm,double *ynorm) 15229b94acceSBarry Smith { 15239b94acceSBarry Smith double norm; 15249b94acceSBarry Smith Scalar cnorm; 1525cddf8d76SBarry Smith VecNorm(y,NORM_2, &norm ); 15269b94acceSBarry Smith if (norm > *delta) { 15279b94acceSBarry Smith norm = *delta/norm; 15289b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 15299b94acceSBarry Smith cnorm = norm; 15309b94acceSBarry Smith VecScale( &cnorm, y ); 15319b94acceSBarry Smith *ynorm = *delta; 15329b94acceSBarry Smith } else { 15339b94acceSBarry Smith *gpnorm = 0.0; 15349b94acceSBarry Smith *ynorm = norm; 15359b94acceSBarry Smith } 15369b94acceSBarry Smith return 0; 15379b94acceSBarry Smith } 15389b94acceSBarry Smith 15395615d1e5SSatish Balay #undef __FUNC__ 15405615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 15419b94acceSBarry Smith /*@ 15429b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 154363a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 15449b94acceSBarry Smith 15459b94acceSBarry Smith Input Parameter: 15469b94acceSBarry Smith . snes - the SNES context 15478ddd3da0SLois Curfman McInnes . x - the solution vector 15489b94acceSBarry Smith 15499b94acceSBarry Smith Output Parameter: 15509b94acceSBarry Smith its - number of iterations until termination 15519b94acceSBarry Smith 15528ddd3da0SLois Curfman McInnes Note: 15538ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 15548ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 15558ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 15568ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 15578ddd3da0SLois Curfman McInnes 15589b94acceSBarry Smith .keywords: SNES, nonlinear, solve 15599b94acceSBarry Smith 156063a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 15619b94acceSBarry Smith @*/ 15628ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 15639b94acceSBarry Smith { 15643c7409f5SSatish Balay int ierr, flg; 1565052efed2SBarry Smith 156677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 156774679c65SBarry Smith PetscValidIntPointer(its); 1568c4fc05e7SBarry Smith if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1569c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 15709b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 1571c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 15729b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 15739b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 15743c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 15756d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 15769b94acceSBarry Smith return 0; 15779b94acceSBarry Smith } 15789b94acceSBarry Smith 15799b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 158037fcc0dbSBarry Smith static NRList *__SNESList = 0; 15819b94acceSBarry Smith 15825615d1e5SSatish Balay #undef __FUNC__ 15835615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 15849b94acceSBarry Smith /*@ 15854b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 15869b94acceSBarry Smith 15879b94acceSBarry Smith Input Parameters: 15889b94acceSBarry Smith . snes - the SNES context 15899b94acceSBarry Smith . method - a known method 15909b94acceSBarry Smith 1591ae12b187SLois Curfman McInnes Options Database Command: 1592ae12b187SLois Curfman McInnes $ -snes_type <method> 1593ae12b187SLois Curfman McInnes $ Use -help for a list of available methods 1594ae12b187SLois Curfman McInnes $ (for instance, ls or tr) 1595ae12b187SLois Curfman McInnes 15969b94acceSBarry Smith Notes: 15979b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 15989b94acceSBarry Smith $ Systems of nonlinear equations: 159940191667SLois Curfman McInnes $ SNES_EQ_LS - Newton's method with line search 160040191667SLois Curfman McInnes $ SNES_EQ_TR - Newton's method with trust region 16019b94acceSBarry Smith $ Unconstrained minimization: 160240191667SLois Curfman McInnes $ SNES_UM_TR - Newton's method with trust region 160340191667SLois Curfman McInnes $ SNES_UM_LS - Newton's method with line search 16049b94acceSBarry Smith 1605ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1606ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1607ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1608ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1609ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1610ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1611ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1612ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1613ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1614ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1615a703fe33SLois Curfman McInnes 1616f536c699SSatish Balay .keywords: SNES, set, method 16179b94acceSBarry Smith @*/ 16184b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 16199b94acceSBarry Smith { 162084cb2905SBarry Smith int ierr; 16219b94acceSBarry Smith int (*r)(SNES); 162277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16237d1a2b2bSBarry Smith if (snes->type == (int) method) return 0; 16247d1a2b2bSBarry Smith 16259b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 162684cb2905SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1627e3372554SBarry Smith if (!__SNESList) {SETERRQ(1,0,"Could not get methods");} 162837fcc0dbSBarry Smith r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1629e3372554SBarry Smith if (!r) {SETERRQ(1,0,"Unknown method");} 16300452661fSBarry Smith if (snes->data) PetscFree(snes->data); 16319b94acceSBarry Smith snes->set_method_called = 1; 16329b94acceSBarry Smith return (*r)(snes); 16339b94acceSBarry Smith } 16349b94acceSBarry Smith 16359b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16365615d1e5SSatish Balay #undef __FUNC__ 1637d4bb536fSBarry Smith #define __FUNC__ "SNESRegister" 16389b94acceSBarry Smith /*@C 16399b94acceSBarry Smith SNESRegister - Adds the method to the nonlinear solver package, given 16404b0e389bSBarry Smith a function pointer and a nonlinear solver name of the type SNESType. 16419b94acceSBarry Smith 16429b94acceSBarry Smith Input Parameters: 16432d872ea7SLois Curfman McInnes . name - either a predefined name such as SNES_EQ_LS, or SNES_NEW 16442d872ea7SLois Curfman McInnes to indicate a new user-defined solver 164540191667SLois Curfman McInnes . sname - corresponding string for name 16469b94acceSBarry Smith . create - routine to create method context 16479b94acceSBarry Smith 164884cb2905SBarry Smith Output Parameter: 164984cb2905SBarry Smith . oname - type associated with this new solver 165084cb2905SBarry Smith 16512d872ea7SLois Curfman McInnes Notes: 16522d872ea7SLois Curfman McInnes Multiple user-defined nonlinear solvers can be added by calling 16532d872ea7SLois Curfman McInnes SNESRegister() with the input parameter "name" set to be SNES_NEW; 16542d872ea7SLois Curfman McInnes each call will return a unique solver type in the output 16552d872ea7SLois Curfman McInnes parameter "oname". 16562d872ea7SLois Curfman McInnes 16579b94acceSBarry Smith .keywords: SNES, nonlinear, register 16589b94acceSBarry Smith 16599b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy() 16609b94acceSBarry Smith @*/ 166184cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES)) 16629b94acceSBarry Smith { 16639b94acceSBarry Smith int ierr; 166484cb2905SBarry Smith static int numberregistered = 0; 166584cb2905SBarry Smith 1666d252947aSBarry Smith if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++); 166784cb2905SBarry Smith 166884cb2905SBarry Smith if (oname) *oname = name; 166937fcc0dbSBarry Smith if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 167084cb2905SBarry Smith NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create ); 16719b94acceSBarry Smith return 0; 16729b94acceSBarry Smith } 1673a847f771SSatish Balay 16749b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16755615d1e5SSatish Balay #undef __FUNC__ 1676d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy" 16779b94acceSBarry Smith /*@C 16789b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 16799b94acceSBarry Smith registered by SNESRegister(). 16809b94acceSBarry Smith 16819b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 16829b94acceSBarry Smith 16839b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 16849b94acceSBarry Smith @*/ 16859b94acceSBarry Smith int SNESRegisterDestroy() 16869b94acceSBarry Smith { 168737fcc0dbSBarry Smith if (__SNESList) { 168837fcc0dbSBarry Smith NRDestroy( __SNESList ); 168937fcc0dbSBarry Smith __SNESList = 0; 16909b94acceSBarry Smith } 169184cb2905SBarry Smith SNESRegisterAllCalled = 0; 16929b94acceSBarry Smith return 0; 16939b94acceSBarry Smith } 16949b94acceSBarry Smith 16955615d1e5SSatish Balay #undef __FUNC__ 1696d4bb536fSBarry Smith #define __FUNC__ "SNESGetTypeFromOptions_Private" 16979b94acceSBarry Smith /* 16984b0e389bSBarry Smith SNESGetTypeFromOptions_Private - Sets the selected method from the 16999b94acceSBarry Smith options database. 17009b94acceSBarry Smith 17019b94acceSBarry Smith Input Parameter: 17029b94acceSBarry Smith . ctx - the SNES context 17039b94acceSBarry Smith 17049b94acceSBarry Smith Output Parameter: 17059b94acceSBarry Smith . method - solver method 17069b94acceSBarry Smith 17079b94acceSBarry Smith Returns: 17089b94acceSBarry Smith Returns 1 if the method is found; 0 otherwise. 17099b94acceSBarry Smith 17109b94acceSBarry Smith Options Database Key: 17114b0e389bSBarry Smith $ -snes_type method 17129b94acceSBarry Smith */ 1713052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 17149b94acceSBarry Smith { 1715052efed2SBarry Smith int ierr; 17169b94acceSBarry Smith char sbuf[50]; 17175baf8537SBarry Smith 1718052efed2SBarry Smith ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1719052efed2SBarry Smith if (*flg) { 1720052efed2SBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 172137fcc0dbSBarry Smith *method = (SNESType)NRFindID( __SNESList, sbuf ); 17229b94acceSBarry Smith } 17239b94acceSBarry Smith return 0; 17249b94acceSBarry Smith } 17259b94acceSBarry Smith 17265615d1e5SSatish Balay #undef __FUNC__ 1727d4bb536fSBarry Smith #define __FUNC__ "SNESGetType" 17289b94acceSBarry Smith /*@C 17299a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 17309b94acceSBarry Smith 17319b94acceSBarry Smith Input Parameter: 17324b0e389bSBarry Smith . snes - nonlinear solver context 17339b94acceSBarry Smith 17349b94acceSBarry Smith Output Parameter: 17359a28b0a6SLois Curfman McInnes . method - SNES method (or use PETSC_NULL) 17369a28b0a6SLois Curfman McInnes . name - name of SNES method (or use PETSC_NULL) 17379b94acceSBarry Smith 17389b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 17399b94acceSBarry Smith @*/ 17404b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name) 17419b94acceSBarry Smith { 174206a9b0adSLois Curfman McInnes int ierr; 174337fcc0dbSBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 17444b0e389bSBarry Smith if (method) *method = (SNESType) snes->type; 174537fcc0dbSBarry Smith if (name) *name = NRFindName( __SNESList, (int) snes->type ); 17469b94acceSBarry Smith return 0; 17479b94acceSBarry Smith } 17489b94acceSBarry Smith 17499b94acceSBarry Smith #include <stdio.h> 17505615d1e5SSatish Balay #undef __FUNC__ 1751d4bb536fSBarry Smith #define __FUNC__ "SNESPrintTypes_Private" 17529b94acceSBarry Smith /* 17534b0e389bSBarry Smith SNESPrintTypes_Private - Prints the SNES methods available from the 17549b94acceSBarry Smith options database. 17559b94acceSBarry Smith 17569b94acceSBarry Smith Input Parameters: 175733455573SSatish Balay . comm - communicator (usually MPI_COMM_WORLD) 17589b94acceSBarry Smith . prefix - prefix (usually "-") 17594b0e389bSBarry Smith . name - the options database name (by default "snes_type") 17609b94acceSBarry Smith */ 176133455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 17629b94acceSBarry Smith { 17639b94acceSBarry Smith FuncList *entry; 176437fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 176537fcc0dbSBarry Smith entry = __SNESList->head; 176677c4ece6SBarry Smith PetscPrintf(comm," %s%s (one of)",prefix,name); 17679b94acceSBarry Smith while (entry) { 176877c4ece6SBarry Smith PetscPrintf(comm," %s",entry->name); 17699b94acceSBarry Smith entry = entry->next; 17709b94acceSBarry Smith } 177177c4ece6SBarry Smith PetscPrintf(comm,"\n"); 17729b94acceSBarry Smith return 0; 17739b94acceSBarry Smith } 17749b94acceSBarry Smith 17755615d1e5SSatish Balay #undef __FUNC__ 1776d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution" 17779b94acceSBarry Smith /*@C 17789b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 17799b94acceSBarry Smith stored. 17809b94acceSBarry Smith 17819b94acceSBarry Smith Input Parameter: 17829b94acceSBarry Smith . snes - the SNES context 17839b94acceSBarry Smith 17849b94acceSBarry Smith Output Parameter: 17859b94acceSBarry Smith . x - the solution 17869b94acceSBarry Smith 17879b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 17889b94acceSBarry Smith 1789a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 17909b94acceSBarry Smith @*/ 17919b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 17929b94acceSBarry Smith { 179377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17949b94acceSBarry Smith *x = snes->vec_sol_always; 17959b94acceSBarry Smith return 0; 17969b94acceSBarry Smith } 17979b94acceSBarry Smith 17985615d1e5SSatish Balay #undef __FUNC__ 1799d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" 18009b94acceSBarry Smith /*@C 18019b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 18029b94acceSBarry Smith stored. 18039b94acceSBarry Smith 18049b94acceSBarry Smith Input Parameter: 18059b94acceSBarry Smith . snes - the SNES context 18069b94acceSBarry Smith 18079b94acceSBarry Smith Output Parameter: 18089b94acceSBarry Smith . x - the solution update 18099b94acceSBarry Smith 18109b94acceSBarry Smith Notes: 18119b94acceSBarry Smith This vector is implementation dependent. 18129b94acceSBarry Smith 18139b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 18149b94acceSBarry Smith 18159b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 18169b94acceSBarry Smith @*/ 18179b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 18189b94acceSBarry Smith { 181977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18209b94acceSBarry Smith *x = snes->vec_sol_update_always; 18219b94acceSBarry Smith return 0; 18229b94acceSBarry Smith } 18239b94acceSBarry Smith 18245615d1e5SSatish Balay #undef __FUNC__ 1825d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction" 18269b94acceSBarry Smith /*@C 18273638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 18289b94acceSBarry Smith 18299b94acceSBarry Smith Input Parameter: 18309b94acceSBarry Smith . snes - the SNES context 18319b94acceSBarry Smith 18329b94acceSBarry Smith Output Parameter: 18333638b69dSLois Curfman McInnes . r - the function 18349b94acceSBarry Smith 18359b94acceSBarry Smith Notes: 18369b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 18379b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 18389b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 18399b94acceSBarry Smith 1840a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 18419b94acceSBarry Smith 184231615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 184331615d04SLois Curfman McInnes SNESGetGradient() 18449b94acceSBarry Smith @*/ 18459b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 18469b94acceSBarry Smith { 184777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1848e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0, 18497c792091SSatish Balay "For SNES_NONLINEAR_EQUATIONS only"); 18509b94acceSBarry Smith *r = snes->vec_func_always; 18519b94acceSBarry Smith return 0; 18529b94acceSBarry Smith } 18539b94acceSBarry Smith 18545615d1e5SSatish Balay #undef __FUNC__ 1855d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient" 18569b94acceSBarry Smith /*@C 18573638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 18589b94acceSBarry Smith 18599b94acceSBarry Smith Input Parameter: 18609b94acceSBarry Smith . snes - the SNES context 18619b94acceSBarry Smith 18629b94acceSBarry Smith Output Parameter: 18639b94acceSBarry Smith . r - the gradient 18649b94acceSBarry Smith 18659b94acceSBarry Smith Notes: 18669b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 18679b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18689b94acceSBarry Smith SNESGetFunction(). 18699b94acceSBarry Smith 18709b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 18719b94acceSBarry Smith 187231615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 18739b94acceSBarry Smith @*/ 18749b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 18759b94acceSBarry Smith { 187677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1877e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 187863c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 18799b94acceSBarry Smith *r = snes->vec_func_always; 18809b94acceSBarry Smith return 0; 18819b94acceSBarry Smith } 18829b94acceSBarry Smith 18835615d1e5SSatish Balay #undef __FUNC__ 1884d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" 18859b94acceSBarry Smith /*@ 18869b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 18879b94acceSBarry Smith unconstrained minimization problems. 18889b94acceSBarry Smith 18899b94acceSBarry Smith Input Parameter: 18909b94acceSBarry Smith . snes - the SNES context 18919b94acceSBarry Smith 18929b94acceSBarry Smith Output Parameter: 18939b94acceSBarry Smith . r - the function 18949b94acceSBarry Smith 18959b94acceSBarry Smith Notes: 18969b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 18979b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18989b94acceSBarry Smith SNESGetFunction(). 18999b94acceSBarry Smith 19009b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 19019b94acceSBarry Smith 190231615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 19039b94acceSBarry Smith @*/ 19049b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 19059b94acceSBarry Smith { 190677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 190774679c65SBarry Smith PetscValidScalarPointer(r); 1908e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 190963c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 19109b94acceSBarry Smith *r = snes->fc; 19119b94acceSBarry Smith return 0; 19129b94acceSBarry Smith } 19139b94acceSBarry Smith 19145615d1e5SSatish Balay #undef __FUNC__ 1915d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" 19163c7409f5SSatish Balay /*@C 19173c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1918d850072dSLois Curfman McInnes SNES options in the database. 19193c7409f5SSatish Balay 19203c7409f5SSatish Balay Input Parameter: 19213c7409f5SSatish Balay . snes - the SNES context 19223c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 19233c7409f5SSatish Balay 1924d850072dSLois Curfman McInnes Notes: 1925a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1926a83b1b31SSatish Balay The first character of all runtime options is AUTOMATICALLY the 1927a83b1b31SSatish Balay hyphen. 1928d850072dSLois Curfman McInnes 19293c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1930a86d99e1SLois Curfman McInnes 1931a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19323c7409f5SSatish Balay @*/ 19333c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 19343c7409f5SSatish Balay { 19353c7409f5SSatish Balay int ierr; 19363c7409f5SSatish Balay 193777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1938639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19393c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19403c7409f5SSatish Balay return 0; 19413c7409f5SSatish Balay } 19423c7409f5SSatish Balay 19435615d1e5SSatish Balay #undef __FUNC__ 1944d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" 19453c7409f5SSatish Balay /*@C 1946f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1947d850072dSLois Curfman McInnes SNES options in the database. 19483c7409f5SSatish Balay 19493c7409f5SSatish Balay Input Parameter: 19503c7409f5SSatish Balay . snes - the SNES context 19513c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 19523c7409f5SSatish Balay 1953d850072dSLois Curfman McInnes Notes: 1954a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1955a83b1b31SSatish Balay The first character of all runtime options is AUTOMATICALLY the 1956a83b1b31SSatish Balay hyphen. 1957d850072dSLois Curfman McInnes 19583c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1959a86d99e1SLois Curfman McInnes 1960a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 19613c7409f5SSatish Balay @*/ 19623c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 19633c7409f5SSatish Balay { 19643c7409f5SSatish Balay int ierr; 19653c7409f5SSatish Balay 196677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1967639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19683c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19693c7409f5SSatish Balay return 0; 19703c7409f5SSatish Balay } 19713c7409f5SSatish Balay 19725615d1e5SSatish Balay #undef __FUNC__ 1973d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" 1974c83590e2SSatish Balay /*@ 19753c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 19763c7409f5SSatish Balay SNES options in the database. 19773c7409f5SSatish Balay 19783c7409f5SSatish Balay Input Parameter: 19793c7409f5SSatish Balay . snes - the SNES context 19803c7409f5SSatish Balay 19813c7409f5SSatish Balay Output Parameter: 19823c7409f5SSatish Balay . prefix - pointer to the prefix string used 19833c7409f5SSatish Balay 19843c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 1985a86d99e1SLois Curfman McInnes 1986a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 19873c7409f5SSatish Balay @*/ 19883c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 19893c7409f5SSatish Balay { 19903c7409f5SSatish Balay int ierr; 19913c7409f5SSatish Balay 199277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1993639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19943c7409f5SSatish Balay return 0; 19953c7409f5SSatish Balay } 19963c7409f5SSatish Balay 19973c7409f5SSatish Balay 19989b94acceSBarry Smith 19999b94acceSBarry Smith 20009b94acceSBarry Smith 2001