19b94acceSBarry Smith #ifndef lint 2*5615d1e5SSatish Balay static char vcid[] = "$Id: snes.c,v 1.104 1997/01/01 03:40:47 bsmith Exp balay $"; 39b94acceSBarry Smith #endif 49b94acceSBarry Smith 59b94acceSBarry Smith #include "draw.h" /*I "draw.h" I*/ 670f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 7f5eb4b81SSatish Balay #include "src/sys/nreg.h" 89b94acceSBarry Smith #include "pinclude/pviewer.h" 99b94acceSBarry Smith #include <math.h> 109b94acceSBarry Smith 11052efed2SBarry Smith extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*); 1233455573SSatish Balay extern int SNESPrintTypes_Private(MPI_Comm,char*,char*); 139b94acceSBarry Smith 14*5615d1e5SSatish Balay #undef __FUNC__ 15*5615d1e5SSatish Balay #define __FUNC__ "SNESView" 169b94acceSBarry Smith /*@ 179b94acceSBarry Smith SNESView - Prints the SNES data structure. 189b94acceSBarry Smith 199b94acceSBarry Smith Input Parameters: 209b94acceSBarry Smith . SNES - the SNES context 219b94acceSBarry Smith . viewer - visualization context 229b94acceSBarry Smith 239b94acceSBarry Smith Options Database Key: 249b94acceSBarry Smith $ -snes_view : calls SNESView() at end of SNESSolve() 259b94acceSBarry Smith 269b94acceSBarry Smith Notes: 279b94acceSBarry Smith The available visualization contexts include 286d4a8577SBarry Smith $ VIEWER_STDOUT_SELF - standard output (default) 296d4a8577SBarry Smith $ VIEWER_STDOUT_WORLD - synchronized standard 309b94acceSBarry Smith $ output where only the first processor opens 319b94acceSBarry Smith $ the file. All other processors send their 329b94acceSBarry Smith $ data to the first processor to print. 339b94acceSBarry Smith 349b94acceSBarry Smith The user can open alternative vistualization contexts with 35dbb450caSBarry Smith $ ViewerFileOpenASCII() - output to a specified file 369b94acceSBarry Smith 379b94acceSBarry Smith .keywords: SNES, view 389b94acceSBarry Smith 39dbb450caSBarry Smith .seealso: ViewerFileOpenASCII() 409b94acceSBarry Smith @*/ 419b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 429b94acceSBarry Smith { 439b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 449b94acceSBarry Smith FILE *fd; 459b94acceSBarry Smith int ierr; 469b94acceSBarry Smith SLES sles; 479b94acceSBarry Smith char *method; 4819bcc07fSBarry Smith ViewerType vtype; 499b94acceSBarry Smith 5074679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5174679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5274679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5374679c65SBarry Smith 5419bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 5519bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 5690ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 5777c4ece6SBarry Smith PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 584b0e389bSBarry Smith SNESGetType(snes,PETSC_NULL,&method); 5977c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," method: %s\n",method); 609b94acceSBarry Smith if (snes->view) (*snes->view)((PetscObject)snes,viewer); 6177c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 629b94acceSBarry Smith " maximum iterations=%d, maximum function evaluations=%d\n", 639b94acceSBarry Smith snes->max_its,snes->max_funcs); 6477c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 659b94acceSBarry Smith " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 669b94acceSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol); 679b94acceSBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 6877c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 699b94acceSBarry Smith if (snes->ksp_ewconv) { 709b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 719b94acceSBarry Smith if (kctx) { 7277c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 739b94acceSBarry Smith " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 749b94acceSBarry Smith kctx->version); 7577c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 769b94acceSBarry Smith " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 779b94acceSBarry Smith kctx->rtol_max,kctx->threshold); 7877c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 799b94acceSBarry Smith kctx->gamma,kctx->alpha,kctx->alpha2); 809b94acceSBarry Smith } 819b94acceSBarry Smith } 82c30f285eSLois Curfman McInnes } else if (vtype == STRING_VIEWER) { 83c30f285eSLois Curfman McInnes SNESGetType(snes,PETSC_NULL,&method); 84c30f285eSLois Curfman McInnes ViewerStringSPrintf(viewer," %-3.3s",method); 8519bcc07fSBarry Smith } 869b94acceSBarry Smith SNESGetSLES(snes,&sles); 879b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 889b94acceSBarry Smith return 0; 899b94acceSBarry Smith } 909b94acceSBarry Smith 91639f9d9dSBarry Smith /* 92639f9d9dSBarry Smith We retain a list of functions that also take SNES command 93639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 94639f9d9dSBarry Smith */ 95639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 96639f9d9dSBarry Smith static int numberofsetfromoptions; 97639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 98639f9d9dSBarry Smith 99*5615d1e5SSatish Balay #undef __FUNC__ 100*5615d1e5SSatish Balay #define __FUNC__ "SNESAddOptionsChecker" 101639f9d9dSBarry Smith /*@ 102639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 103639f9d9dSBarry Smith 104639f9d9dSBarry Smith Input Parameter: 105639f9d9dSBarry Smith . snescheck - function that checks for options 106639f9d9dSBarry Smith 107639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 108639f9d9dSBarry Smith @*/ 109639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 110639f9d9dSBarry Smith { 111639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 112e3372554SBarry Smith SETERRQ(1,0,"Too many options checkers, only 5 allowed"); 113639f9d9dSBarry Smith } 114639f9d9dSBarry Smith 115639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 116639f9d9dSBarry Smith return 0; 117639f9d9dSBarry Smith } 118639f9d9dSBarry Smith 119*5615d1e5SSatish Balay #undef __FUNC__ 120*5615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1219b94acceSBarry Smith /*@ 122682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1239b94acceSBarry Smith 1249b94acceSBarry Smith Input Parameter: 1259b94acceSBarry Smith . snes - the SNES context 1269b94acceSBarry Smith 1279b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1289b94acceSBarry Smith 129a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1309b94acceSBarry Smith @*/ 1319b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1329b94acceSBarry Smith { 1334b0e389bSBarry Smith SNESType method; 1349b94acceSBarry Smith double tmp; 1359b94acceSBarry Smith SLES sles; 13681f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1373c7409f5SSatish Balay int version = PETSC_DEFAULT; 1389b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1399b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1409b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1419b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1429b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1439b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1449b94acceSBarry Smith 14581f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 14681f57ec7SBarry Smith 14781f57ec7SBarry Smith 14877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 149e3372554SBarry Smith if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp"); 150052efed2SBarry Smith ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 151052efed2SBarry Smith if (flg) { 152052efed2SBarry Smith ierr = SNESSetType(snes,method); CHKERRQ(ierr); 1539b94acceSBarry Smith } 1545334005bSBarry Smith else if (!snes->set_method_called) { 1555334005bSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 15640191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 1575334005bSBarry Smith } 1585334005bSBarry Smith else { 15940191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 1605334005bSBarry Smith } 1615334005bSBarry Smith } 1625334005bSBarry Smith 1633c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1643c7409f5SSatish Balay if (flg) { SNESPrintHelp(snes); } 1653c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 166d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);} 1673c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 168d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);} 1693c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 170d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); } 1713c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 1723c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 173d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 174d7a720efSLois Curfman McInnes if (flg) { SNESSetTrustRegionTolerance(snes,tmp); } 175d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 176d7a720efSLois Curfman McInnes if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); } 1773c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 1783c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 179b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 180b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 181b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 182b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 183b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 184b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 185b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 1869b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 1879b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 1885bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 1895bbad29bSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);} 1903c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 191639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 1923c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 193d31a3109SLois Curfman McInnes if (flg) { 194639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 195d31a3109SLois Curfman McInnes } 19681f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 1973c7409f5SSatish Balay if (flg) { 19817699dbbSLois Curfman McInnes int rank = 0; 199d7e8b826SBarry Smith DrawLG lg; 20017699dbbSLois Curfman McInnes MPI_Initialized(&rank); 20117699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 20217699dbbSLois Curfman McInnes if (!rank) { 20381f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2049b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 205f8c826e1SBarry Smith snes->xmonitor = lg; 2069b94acceSBarry Smith PLogObjectParent(snes,lg); 2079b94acceSBarry Smith } 2089b94acceSBarry Smith } 2093c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2103c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2119b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2129b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 21394a424c1SBarry Smith PLogInfo(snes, 21431615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 21531615d04SLois Curfman McInnes } 21631615d04SLois Curfman McInnes else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 21731615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 21831615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 21994a424c1SBarry Smith PLogInfo(snes, 22031615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2219b94acceSBarry Smith } 222639f9d9dSBarry Smith 223639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 224639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 225639f9d9dSBarry Smith } 226639f9d9dSBarry Smith 2279b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2289b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 2299b94acceSBarry Smith if (!snes->setfromoptions) return 0; 2309b94acceSBarry Smith return (*snes->setfromoptions)(snes); 2319b94acceSBarry Smith } 2329b94acceSBarry Smith 233*5615d1e5SSatish Balay #undef __FUNC__ 234*5615d1e5SSatish Balay #define __FUNC__ "SNESPrintHelp" 2359b94acceSBarry Smith /*@ 2369b94acceSBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2379b94acceSBarry Smith 2389b94acceSBarry Smith Input Parameter: 2399b94acceSBarry Smith . snes - the SNES context 2409b94acceSBarry Smith 241a703fe33SLois Curfman McInnes Options Database Keys: 242a703fe33SLois Curfman McInnes $ -help, -h 243a703fe33SLois Curfman McInnes 2449b94acceSBarry Smith .keywords: SNES, nonlinear, help 2459b94acceSBarry Smith 246682d7d0cSBarry Smith .seealso: SNESSetFromOptions() 2479b94acceSBarry Smith @*/ 2489b94acceSBarry Smith int SNESPrintHelp(SNES snes) 2499b94acceSBarry Smith { 2506daaf66cSBarry Smith char p[64]; 2516daaf66cSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2526daaf66cSBarry Smith 25377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2546daaf66cSBarry Smith 2556daaf66cSBarry Smith PetscStrcpy(p,"-"); 2566daaf66cSBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2576daaf66cSBarry Smith 2586daaf66cSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2596daaf66cSBarry Smith 260d31a3109SLois Curfman McInnes PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 26133455573SSatish Balay SNESPrintTypes_Private(snes->comm,p,"snes_type"); 26277c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 263b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 264b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 265b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 266b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 267b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 268b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 2695bbad29bSBarry Smith PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 2705bbad29bSBarry Smith PetscPrintf(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 271d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor\n",p); 272d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 273b18e04deSLois Curfman McInnes if (snes->type == SNES_NONLINEAR_EQUATIONS) { 27477c4ece6SBarry Smith PetscPrintf(snes->comm, 275d31a3109SLois Curfman McInnes " Options for solving systems of nonlinear equations only:\n"); 27677c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 27777c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 278ef1dfb11SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 2791650c7b0SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err: relative error in function evaluation for matrix-free Jacobian\n",p); 2801650c7b0SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Jacobian\n",p); 28177c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 28277c4ece6SBarry Smith PetscPrintf(snes->comm, 283b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 28477c4ece6SBarry Smith PetscPrintf(snes->comm, 285b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 28677c4ece6SBarry Smith PetscPrintf(snes->comm, 287b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 28877c4ece6SBarry Smith PetscPrintf(snes->comm, 289b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 29077c4ece6SBarry Smith PetscPrintf(snes->comm, 291b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 29277c4ece6SBarry Smith PetscPrintf(snes->comm, 293b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 29477c4ece6SBarry Smith PetscPrintf(snes->comm, 295b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 296b18e04deSLois Curfman McInnes } 297b18e04deSLois Curfman McInnes else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 29877c4ece6SBarry Smith PetscPrintf(snes->comm, 299d31a3109SLois Curfman McInnes " Options for solving unconstrained minimization problems only:\n"); 300b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 30177c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 30277c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 303d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_err: relative error in gradient evaluation for matrix-free Hessian\n",p); 304d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_umin: minimum iterate parameter for matrix-free Hessian\n",p); 305b18e04deSLois Curfman McInnes } 3064537a946SLois Curfman McInnes PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 30777c4ece6SBarry Smith PetscPrintf(snes->comm,"a particular method\n"); 3086daaf66cSBarry Smith if (snes->printhelp) (*snes->printhelp)(snes,p); 3099b94acceSBarry Smith return 0; 3109b94acceSBarry Smith } 311a847f771SSatish Balay 312*5615d1e5SSatish Balay #undef __FUNC__ 313*5615d1e5SSatish Balay #define __FUNC__ "SNESSetApplicationContext" 3149b94acceSBarry Smith /*@ 3159b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 3169b94acceSBarry Smith the nonlinear solvers. 3179b94acceSBarry Smith 3189b94acceSBarry Smith Input Parameters: 3199b94acceSBarry Smith . snes - the SNES context 3209b94acceSBarry Smith . usrP - optional user context 3219b94acceSBarry Smith 3229b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3239b94acceSBarry Smith 3249b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3259b94acceSBarry Smith @*/ 3269b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3279b94acceSBarry Smith { 32877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3299b94acceSBarry Smith snes->user = usrP; 3309b94acceSBarry Smith return 0; 3319b94acceSBarry Smith } 33274679c65SBarry Smith 333*5615d1e5SSatish Balay #undef __FUNC__ 334*5615d1e5SSatish Balay #define __FUNC__ "SNESGetApplicationContext" 3359b94acceSBarry Smith /*@C 3369b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3379b94acceSBarry Smith nonlinear solvers. 3389b94acceSBarry Smith 3399b94acceSBarry Smith Input Parameter: 3409b94acceSBarry Smith . snes - SNES context 3419b94acceSBarry Smith 3429b94acceSBarry Smith Output Parameter: 3439b94acceSBarry Smith . usrP - user context 3449b94acceSBarry Smith 3459b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3469b94acceSBarry Smith 3479b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3489b94acceSBarry Smith @*/ 3499b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3509b94acceSBarry Smith { 35177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3529b94acceSBarry Smith *usrP = snes->user; 3539b94acceSBarry Smith return 0; 3549b94acceSBarry Smith } 35574679c65SBarry Smith 356*5615d1e5SSatish Balay #undef __FUNC__ 357*5615d1e5SSatish Balay #define __FUNC__ "SNESGetIterationNumber" 3589b94acceSBarry Smith /*@ 3599b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3609b94acceSBarry Smith nonlinear solver. 3619b94acceSBarry Smith 3629b94acceSBarry Smith Input Parameter: 3639b94acceSBarry Smith . snes - SNES context 3649b94acceSBarry Smith 3659b94acceSBarry Smith Output Parameter: 3669b94acceSBarry Smith . iter - iteration number 3679b94acceSBarry Smith 3689b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3699b94acceSBarry Smith @*/ 3709b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3719b94acceSBarry Smith { 37277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 37374679c65SBarry Smith PetscValidIntPointer(iter); 3749b94acceSBarry Smith *iter = snes->iter; 3759b94acceSBarry Smith return 0; 3769b94acceSBarry Smith } 37774679c65SBarry Smith 378*5615d1e5SSatish Balay #undef __FUNC__ 379*5615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 3809b94acceSBarry Smith /*@ 3819b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3829b94acceSBarry Smith with SNESSSetFunction(). 3839b94acceSBarry Smith 3849b94acceSBarry Smith Input Parameter: 3859b94acceSBarry Smith . snes - SNES context 3869b94acceSBarry Smith 3879b94acceSBarry Smith Output Parameter: 3889b94acceSBarry Smith . fnorm - 2-norm of function 3899b94acceSBarry Smith 3909b94acceSBarry Smith Note: 3919b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 392a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 393a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 3949b94acceSBarry Smith 3959b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 396a86d99e1SLois Curfman McInnes 397a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 3989b94acceSBarry Smith @*/ 3999b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 4009b94acceSBarry Smith { 40177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 40274679c65SBarry Smith PetscValidScalarPointer(fnorm); 40374679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 404e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 40574679c65SBarry Smith } 4069b94acceSBarry Smith *fnorm = snes->norm; 4079b94acceSBarry Smith return 0; 4089b94acceSBarry Smith } 40974679c65SBarry Smith 410*5615d1e5SSatish Balay #undef __FUNC__ 411*5615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 4129b94acceSBarry Smith /*@ 4139b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 4149b94acceSBarry Smith with SNESSSetGradient(). 4159b94acceSBarry Smith 4169b94acceSBarry Smith Input Parameter: 4179b94acceSBarry Smith . snes - SNES context 4189b94acceSBarry Smith 4199b94acceSBarry Smith Output Parameter: 4209b94acceSBarry Smith . fnorm - 2-norm of gradient 4219b94acceSBarry Smith 4229b94acceSBarry Smith Note: 4239b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 424a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 425a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4269b94acceSBarry Smith 4279b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 428a86d99e1SLois Curfman McInnes 429a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4309b94acceSBarry Smith @*/ 4319b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4329b94acceSBarry Smith { 43377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 43474679c65SBarry Smith PetscValidScalarPointer(gnorm); 43574679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 436e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 43774679c65SBarry Smith } 4389b94acceSBarry Smith *gnorm = snes->norm; 4399b94acceSBarry Smith return 0; 4409b94acceSBarry Smith } 44174679c65SBarry Smith 442*5615d1e5SSatish Balay #undef __FUNC__ 443*5615d1e5SSatish Balay #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 4449b94acceSBarry Smith /*@ 4459b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4469b94acceSBarry Smith attempted by the nonlinear solver. 4479b94acceSBarry Smith 4489b94acceSBarry Smith Input Parameter: 4499b94acceSBarry Smith . snes - SNES context 4509b94acceSBarry Smith 4519b94acceSBarry Smith Output Parameter: 4529b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4539b94acceSBarry Smith 4549b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4559b94acceSBarry Smith @*/ 4569b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4579b94acceSBarry Smith { 45877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 45974679c65SBarry Smith PetscValidIntPointer(nfails); 4609b94acceSBarry Smith *nfails = snes->nfailures; 4619b94acceSBarry Smith return 0; 4629b94acceSBarry Smith } 463a847f771SSatish Balay 464*5615d1e5SSatish Balay #undef __FUNC__ 465*5615d1e5SSatish Balay #define __FUNC__ "SNESGetSLES" 4669b94acceSBarry Smith /*@C 4679b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 4689b94acceSBarry Smith 4699b94acceSBarry Smith Input Parameter: 4709b94acceSBarry Smith . snes - the SNES context 4719b94acceSBarry Smith 4729b94acceSBarry Smith Output Parameter: 4739b94acceSBarry Smith . sles - the SLES context 4749b94acceSBarry Smith 4759b94acceSBarry Smith Notes: 4769b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 4779b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 4789b94acceSBarry Smith KSP and PC contexts as well. 4799b94acceSBarry Smith 4809b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 4819b94acceSBarry Smith 4829b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 4839b94acceSBarry Smith @*/ 4849b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 4859b94acceSBarry Smith { 48677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 4879b94acceSBarry Smith *sles = snes->sles; 4889b94acceSBarry Smith return 0; 4899b94acceSBarry Smith } 4909b94acceSBarry Smith /* -----------------------------------------------------------*/ 491*5615d1e5SSatish Balay #undef __FUNC__ 492*5615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 4939b94acceSBarry Smith /*@C 4949b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 4959b94acceSBarry Smith 4969b94acceSBarry Smith Input Parameter: 4979b94acceSBarry Smith . comm - MPI communicator 4989b94acceSBarry Smith . type - type of method, one of 4999b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS 5009b94acceSBarry Smith $ (for systems of nonlinear equations) 5019b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION 5029b94acceSBarry Smith $ (for unconstrained minimization) 5039b94acceSBarry Smith 5049b94acceSBarry Smith Output Parameter: 5059b94acceSBarry Smith . outsnes - the new SNES context 5069b94acceSBarry Smith 507c1f60f51SBarry Smith Options Database Key: 50819bd6747SLois Curfman McInnes $ -snes_mf - use default matrix-free Jacobian-vector products, 50919bd6747SLois Curfman McInnes $ and no preconditioning matrix 51019bd6747SLois Curfman McInnes $ -snes_mf_operator - use default matrix-free Jacobian-vector 51119bd6747SLois Curfman McInnes $ products, and a user-provided preconditioning matrix 51219bd6747SLois Curfman McInnes $ as set by SNESSetJacobian() 51319bd6747SLois Curfman McInnes $ -snes_fd - use (slow!) finite differences to compute Jacobian 514c1f60f51SBarry Smith 5159b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5169b94acceSBarry Smith 51763a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 5189b94acceSBarry Smith @*/ 5194b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 5209b94acceSBarry Smith { 5219b94acceSBarry Smith int ierr; 5229b94acceSBarry Smith SNES snes; 5239b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 52437fcc0dbSBarry Smith 5259b94acceSBarry Smith *outsnes = 0; 5262a0acf01SLois Curfman McInnes if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 527e3372554SBarry Smith SETERRQ(1,0,"incorrect method type"); 5280452661fSBarry Smith PetscHeaderCreate(snes,_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm); 5299b94acceSBarry Smith PLogObjectCreate(snes); 5309b94acceSBarry Smith snes->max_its = 50; 5319b94acceSBarry Smith snes->max_funcs = 1000; 5329b94acceSBarry Smith snes->norm = 0.0; 5335a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 534b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 535b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5369b94acceSBarry Smith snes->atol = 1.e-10; 5375a2d0531SBarry Smith } 538b4874afaSBarry Smith else { 539b4874afaSBarry Smith snes->rtol = 1.e-8; 540b4874afaSBarry Smith snes->ttol = 0.0; 541b4874afaSBarry Smith snes->atol = 1.e-50; 542b4874afaSBarry Smith } 5439b94acceSBarry Smith snes->xtol = 1.e-8; 544b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5459b94acceSBarry Smith snes->nfuncs = 0; 5469b94acceSBarry Smith snes->nfailures = 0; 547639f9d9dSBarry Smith snes->numbermonitors = 0; 5489b94acceSBarry Smith snes->data = 0; 5499b94acceSBarry Smith snes->view = 0; 5509b94acceSBarry Smith snes->computeumfunction = 0; 5519b94acceSBarry Smith snes->umfunP = 0; 5529b94acceSBarry Smith snes->fc = 0; 5539b94acceSBarry Smith snes->deltatol = 1.e-12; 5549b94acceSBarry Smith snes->fmin = -1.e30; 5559b94acceSBarry Smith snes->method_class = type; 5569b94acceSBarry Smith snes->set_method_called = 0; 557a703fe33SLois Curfman McInnes snes->setup_called = 0; 5589b94acceSBarry Smith snes->ksp_ewconv = 0; 5597d1a2b2bSBarry Smith snes->type = -1; 5606f24a144SLois Curfman McInnes snes->xmonitor = 0; 5616f24a144SLois Curfman McInnes snes->vwork = 0; 5626f24a144SLois Curfman McInnes snes->nwork = 0; 5639b94acceSBarry Smith 5649b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 5650452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 5669b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 5679b94acceSBarry Smith kctx->version = 2; 5689b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 5699b94acceSBarry Smith this was too large for some test cases */ 5709b94acceSBarry Smith kctx->rtol_last = 0; 5719b94acceSBarry Smith kctx->rtol_max = .9; 5729b94acceSBarry Smith kctx->gamma = 1.0; 5739b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 5749b94acceSBarry Smith kctx->alpha = kctx->alpha2; 5759b94acceSBarry Smith kctx->threshold = .1; 5769b94acceSBarry Smith kctx->lresid_last = 0; 5779b94acceSBarry Smith kctx->norm_last = 0; 5789b94acceSBarry Smith 5799b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 5809b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 5815334005bSBarry Smith 5829b94acceSBarry Smith *outsnes = snes; 5839b94acceSBarry Smith return 0; 5849b94acceSBarry Smith } 5859b94acceSBarry Smith 5869b94acceSBarry Smith /* --------------------------------------------------------------- */ 587*5615d1e5SSatish Balay #undef __FUNC__ 588*5615d1e5SSatish Balay #define __FUNC__ "SNESSetFunction" 5899b94acceSBarry Smith /*@C 5909b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 5919b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 5929b94acceSBarry Smith equations. 5939b94acceSBarry Smith 5949b94acceSBarry Smith Input Parameters: 5959b94acceSBarry Smith . snes - the SNES context 5969b94acceSBarry Smith . func - function evaluation routine 5979b94acceSBarry Smith . r - vector to store function value 5982cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 5992cd2dfdcSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6009b94acceSBarry Smith 6019b94acceSBarry Smith Calling sequence of func: 602313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Vec f,void *ctx); 6039b94acceSBarry Smith 6049b94acceSBarry Smith . x - input vector 605313e4042SLois Curfman McInnes . f - function vector 6062cd2dfdcSLois Curfman McInnes . ctx - optional user-defined function context 6079b94acceSBarry Smith 6089b94acceSBarry Smith Notes: 6099b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6109b94acceSBarry Smith $ f'(x) x = -f(x), 6119b94acceSBarry Smith $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 6129b94acceSBarry Smith 6139b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6149b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6159b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 6169b94acceSBarry Smith 6179b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6189b94acceSBarry Smith 619a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6209b94acceSBarry Smith @*/ 6215334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 6229b94acceSBarry Smith { 62377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 624e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6259b94acceSBarry Smith snes->computefunction = func; 6269b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6279b94acceSBarry Smith snes->funP = ctx; 6289b94acceSBarry Smith return 0; 6299b94acceSBarry Smith } 6309b94acceSBarry Smith 631*5615d1e5SSatish Balay #undef __FUNC__ 632*5615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 6339b94acceSBarry Smith /*@ 6349b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6359b94acceSBarry Smith SNESSetFunction(). 6369b94acceSBarry Smith 6379b94acceSBarry Smith Input Parameters: 6389b94acceSBarry Smith . snes - the SNES context 6399b94acceSBarry Smith . x - input vector 6409b94acceSBarry Smith 6419b94acceSBarry Smith Output Parameter: 6423638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6439b94acceSBarry Smith 6441bffabb2SLois Curfman McInnes Notes: 6459b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6469b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6479b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6489b94acceSBarry Smith 6499b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6509b94acceSBarry Smith 651a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6529b94acceSBarry Smith @*/ 6539b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 6549b94acceSBarry Smith { 6559b94acceSBarry Smith int ierr; 6569b94acceSBarry Smith 65774679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 658e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6599b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 6609b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 6619b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 6629b94acceSBarry Smith return 0; 6639b94acceSBarry Smith } 6649b94acceSBarry Smith 665*5615d1e5SSatish Balay #undef __FUNC__ 666*5615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunction" 6679b94acceSBarry Smith /*@C 6689b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 6699b94acceSBarry Smith unconstrained minimization. 6709b94acceSBarry Smith 6719b94acceSBarry Smith Input Parameters: 6729b94acceSBarry Smith . snes - the SNES context 6739b94acceSBarry Smith . func - function evaluation routine 674044dda88SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 675044dda88SLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6769b94acceSBarry Smith 6779b94acceSBarry Smith Calling sequence of func: 6789b94acceSBarry Smith . func (SNES snes,Vec x,double *f,void *ctx); 6799b94acceSBarry Smith 6809b94acceSBarry Smith . x - input vector 6819b94acceSBarry Smith . f - function 682044dda88SLois Curfman McInnes . ctx - [optional] user-defined function context 6839b94acceSBarry Smith 6849b94acceSBarry Smith Notes: 6859b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 6869b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 6879b94acceSBarry Smith SNESSetFunction(). 6889b94acceSBarry Smith 6899b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 6909b94acceSBarry Smith 691a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 692a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 6939b94acceSBarry Smith @*/ 6949b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 6959b94acceSBarry Smith void *ctx) 6969b94acceSBarry Smith { 69777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 698e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 6999b94acceSBarry Smith snes->computeumfunction = func; 7009b94acceSBarry Smith snes->umfunP = ctx; 7019b94acceSBarry Smith return 0; 7029b94acceSBarry Smith } 7039b94acceSBarry Smith 704*5615d1e5SSatish Balay #undef __FUNC__ 705*5615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 7069b94acceSBarry Smith /*@ 7079b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 7089b94acceSBarry Smith set with SNESSetMinimizationFunction(). 7099b94acceSBarry Smith 7109b94acceSBarry Smith Input Parameters: 7119b94acceSBarry Smith . snes - the SNES context 7129b94acceSBarry Smith . x - input vector 7139b94acceSBarry Smith 7149b94acceSBarry Smith Output Parameter: 7159b94acceSBarry Smith . y - function value 7169b94acceSBarry Smith 7179b94acceSBarry Smith Notes: 7189b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 7199b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7209b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 721a86d99e1SLois Curfman McInnes 722a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 723a86d99e1SLois Curfman McInnes 724a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 725a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 7269b94acceSBarry Smith @*/ 7279b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 7289b94acceSBarry Smith { 7299b94acceSBarry Smith int ierr; 730e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7319b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 7329b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 7339b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7349b94acceSBarry Smith return 0; 7359b94acceSBarry Smith } 7369b94acceSBarry Smith 737*5615d1e5SSatish Balay #undef __FUNC__ 738*5615d1e5SSatish Balay #define __FUNC__ "SNESSetGradient" 7399b94acceSBarry Smith /*@C 7409b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 7419b94acceSBarry Smith vector for use by the SNES routines. 7429b94acceSBarry Smith 7439b94acceSBarry Smith Input Parameters: 7449b94acceSBarry Smith . snes - the SNES context 7459b94acceSBarry Smith . func - function evaluation routine 746044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 747044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 7489b94acceSBarry Smith . r - vector to store gradient value 7499b94acceSBarry Smith 7509b94acceSBarry Smith Calling sequence of func: 7519b94acceSBarry Smith . func (SNES, Vec x, Vec g, void *ctx); 7529b94acceSBarry Smith 7539b94acceSBarry Smith . x - input vector 7549b94acceSBarry Smith . g - gradient vector 755044dda88SLois Curfman McInnes . ctx - optional user-defined gradient context 7569b94acceSBarry Smith 7579b94acceSBarry Smith Notes: 7589b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7599b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7609b94acceSBarry Smith SNESSetFunction(). 7619b94acceSBarry Smith 7629b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 7639b94acceSBarry Smith 764a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 765a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 7669b94acceSBarry Smith @*/ 76774679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 7689b94acceSBarry Smith { 76977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 770e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 7719b94acceSBarry Smith snes->computefunction = func; 7729b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7739b94acceSBarry Smith snes->funP = ctx; 7749b94acceSBarry Smith return 0; 7759b94acceSBarry Smith } 7769b94acceSBarry Smith 777*5615d1e5SSatish Balay #undef __FUNC__ 778*5615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 7799b94acceSBarry Smith /*@ 780a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 781a86d99e1SLois Curfman McInnes SNESSetGradient(). 7829b94acceSBarry Smith 7839b94acceSBarry Smith Input Parameters: 7849b94acceSBarry Smith . snes - the SNES context 7859b94acceSBarry Smith . x - input vector 7869b94acceSBarry Smith 7879b94acceSBarry Smith Output Parameter: 7889b94acceSBarry Smith . y - gradient vector 7899b94acceSBarry Smith 7909b94acceSBarry Smith Notes: 7919b94acceSBarry Smith SNESComputeGradient() is valid only for 7929b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7939b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 794a86d99e1SLois Curfman McInnes 795a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 796a86d99e1SLois Curfman McInnes 797a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 798a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 7999b94acceSBarry Smith @*/ 8009b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 8019b94acceSBarry Smith { 8029b94acceSBarry Smith int ierr; 80374679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 804e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8059b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 8069b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 8079b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 8089b94acceSBarry Smith return 0; 8099b94acceSBarry Smith } 8109b94acceSBarry Smith 811*5615d1e5SSatish Balay #undef __FUNC__ 812*5615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 81362fef451SLois Curfman McInnes /*@ 81462fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 81562fef451SLois Curfman McInnes set with SNESSetJacobian(). 81662fef451SLois Curfman McInnes 81762fef451SLois Curfman McInnes Input Parameters: 81862fef451SLois Curfman McInnes . snes - the SNES context 81962fef451SLois Curfman McInnes . x - input vector 82062fef451SLois Curfman McInnes 82162fef451SLois Curfman McInnes Output Parameters: 82262fef451SLois Curfman McInnes . A - Jacobian matrix 82362fef451SLois Curfman McInnes . B - optional preconditioning matrix 82462fef451SLois Curfman McInnes . flag - flag indicating matrix structure 82562fef451SLois Curfman McInnes 82662fef451SLois Curfman McInnes Notes: 82762fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 82862fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 82962fef451SLois Curfman McInnes 830dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 831dc5a77f8SLois Curfman McInnes flag parameter. 83262fef451SLois Curfman McInnes 83362fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 83462fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 83562fef451SLois Curfman McInnes methods is SNESComputeJacobian(). 83662fef451SLois Curfman McInnes 83762fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 83862fef451SLois Curfman McInnes 83962fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 84062fef451SLois Curfman McInnes @*/ 8419b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8429b94acceSBarry Smith { 8439b94acceSBarry Smith int ierr; 84474679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 845e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 8469b94acceSBarry Smith if (!snes->computejacobian) return 0; 8479b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 848c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 8499b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 8509b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 8516d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 85277c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 85377c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 8549b94acceSBarry Smith return 0; 8559b94acceSBarry Smith } 8569b94acceSBarry Smith 857*5615d1e5SSatish Balay #undef __FUNC__ 858*5615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 85962fef451SLois Curfman McInnes /*@ 86062fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 86162fef451SLois Curfman McInnes set with SNESSetHessian(). 86262fef451SLois Curfman McInnes 86362fef451SLois Curfman McInnes Input Parameters: 86462fef451SLois Curfman McInnes . snes - the SNES context 86562fef451SLois Curfman McInnes . x - input vector 86662fef451SLois Curfman McInnes 86762fef451SLois Curfman McInnes Output Parameters: 86862fef451SLois Curfman McInnes . A - Hessian matrix 86962fef451SLois Curfman McInnes . B - optional preconditioning matrix 87062fef451SLois Curfman McInnes . flag - flag indicating matrix structure 87162fef451SLois Curfman McInnes 87262fef451SLois Curfman McInnes Notes: 87362fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 87462fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 87562fef451SLois Curfman McInnes 876dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 877dc5a77f8SLois Curfman McInnes flag parameter. 87862fef451SLois Curfman McInnes 87962fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 88062fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 88162fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 88262fef451SLois Curfman McInnes 88362fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 88462fef451SLois Curfman McInnes 885a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 886a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 88762fef451SLois Curfman McInnes @*/ 88862fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 8899b94acceSBarry Smith { 8909b94acceSBarry Smith int ierr; 89174679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 892e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8939b94acceSBarry Smith if (!snes->computejacobian) return 0; 89462fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 8950f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 89662fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 89762fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 8980f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 89977c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 90077c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9019b94acceSBarry Smith return 0; 9029b94acceSBarry Smith } 9039b94acceSBarry Smith 904*5615d1e5SSatish Balay #undef __FUNC__ 905*5615d1e5SSatish Balay #define __FUNC__ "SNESSetJacobian" 9069b94acceSBarry Smith /*@C 9079b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 908044dda88SLois Curfman McInnes location to store the matrix. 9099b94acceSBarry Smith 9109b94acceSBarry Smith Input Parameters: 9119b94acceSBarry Smith . snes - the SNES context 9129b94acceSBarry Smith . A - Jacobian matrix 9139b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 9149b94acceSBarry Smith . func - Jacobian evaluation routine 9152cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 9162cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 9179b94acceSBarry Smith 9189b94acceSBarry Smith Calling sequence of func: 919313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9209b94acceSBarry Smith 9219b94acceSBarry Smith . x - input vector 9229b94acceSBarry Smith . A - Jacobian matrix 9239b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 924ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 925ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 9262cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Jacobian context 9279b94acceSBarry Smith 9289b94acceSBarry Smith Notes: 929dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 9302cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 931ac21db08SLois Curfman McInnes 932ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9339b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9349b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9359b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9369b94acceSBarry Smith throughout the global iterations. 9379b94acceSBarry Smith 9389b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9399b94acceSBarry Smith 940ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 9419b94acceSBarry Smith @*/ 9429b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 9439b94acceSBarry Smith MatStructure*,void*),void *ctx) 9449b94acceSBarry Smith { 94577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 94674679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 947e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 9489b94acceSBarry Smith snes->computejacobian = func; 9499b94acceSBarry Smith snes->jacP = ctx; 9509b94acceSBarry Smith snes->jacobian = A; 9519b94acceSBarry Smith snes->jacobian_pre = B; 9529b94acceSBarry Smith return 0; 9539b94acceSBarry Smith } 95462fef451SLois Curfman McInnes 955*5615d1e5SSatish Balay #undef __FUNC__ 956*5615d1e5SSatish Balay #define __FUNC__ "SNESGetJacobian" 957b4fd4287SBarry Smith /*@ 958b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 959b4fd4287SBarry Smith provided context for evaluating the Jacobian. 960b4fd4287SBarry Smith 961b4fd4287SBarry Smith Input Parameter: 962b4fd4287SBarry Smith . snes - the nonlinear solver context 963b4fd4287SBarry Smith 964b4fd4287SBarry Smith Output Parameters: 965b4fd4287SBarry Smith . A - location to stash Jacobian matrix (or PETSC_NULL) 966b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 967b4fd4287SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 968b4fd4287SBarry Smith 969b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 970b4fd4287SBarry Smith @*/ 971b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 972b4fd4287SBarry Smith { 97374679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 974e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 975b4fd4287SBarry Smith if (A) *A = snes->jacobian; 976b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 977b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 978b4fd4287SBarry Smith return 0; 979b4fd4287SBarry Smith } 980b4fd4287SBarry Smith 981*5615d1e5SSatish Balay #undef __FUNC__ 982*5615d1e5SSatish Balay #define __FUNC__ "SNESSetHessian" 9839b94acceSBarry Smith /*@C 9849b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 985044dda88SLois Curfman McInnes location to store the matrix. 9869b94acceSBarry Smith 9879b94acceSBarry Smith Input Parameters: 9889b94acceSBarry Smith . snes - the SNES context 9899b94acceSBarry Smith . A - Hessian matrix 9909b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 9919b94acceSBarry Smith . func - Jacobian evaluation routine 992313e4042SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 993313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 9949b94acceSBarry Smith 9959b94acceSBarry Smith Calling sequence of func: 996313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9979b94acceSBarry Smith 9989b94acceSBarry Smith . x - input vector 9999b94acceSBarry Smith . A - Hessian matrix 10009b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1001ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1002ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 10032cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Hessian context 10049b94acceSBarry Smith 10059b94acceSBarry Smith Notes: 1006dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10072cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1008ac21db08SLois Curfman McInnes 10099b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 10109b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 10119b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10129b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10139b94acceSBarry Smith throughout the global iterations. 10149b94acceSBarry Smith 10159b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 10169b94acceSBarry Smith 1017ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 10189b94acceSBarry Smith @*/ 10199b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 10209b94acceSBarry Smith MatStructure*,void*),void *ctx) 10219b94acceSBarry Smith { 102277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 102374679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1024e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 10259b94acceSBarry Smith snes->computejacobian = func; 10269b94acceSBarry Smith snes->jacP = ctx; 10279b94acceSBarry Smith snes->jacobian = A; 10289b94acceSBarry Smith snes->jacobian_pre = B; 10299b94acceSBarry Smith return 0; 10309b94acceSBarry Smith } 10319b94acceSBarry Smith 1032*5615d1e5SSatish Balay #undef __FUNC__ 1033*5615d1e5SSatish Balay #define __FUNC__ "SNESGetHessian" 103462fef451SLois Curfman McInnes /*@ 103562fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 103662fef451SLois Curfman McInnes provided context for evaluating the Hessian. 103762fef451SLois Curfman McInnes 103862fef451SLois Curfman McInnes Input Parameter: 103962fef451SLois Curfman McInnes . snes - the nonlinear solver context 104062fef451SLois Curfman McInnes 104162fef451SLois Curfman McInnes Output Parameters: 104262fef451SLois Curfman McInnes . A - location to stash Hessian matrix (or PETSC_NULL) 104362fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 104462fef451SLois Curfman McInnes . ctx - location to stash Hessian ctx (or PETSC_NULL) 104562fef451SLois Curfman McInnes 104662fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 104762fef451SLois Curfman McInnes @*/ 104862fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 104962fef451SLois Curfman McInnes { 105074679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1051e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 105262fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 105362fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 105462fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 105562fef451SLois Curfman McInnes return 0; 105662fef451SLois Curfman McInnes } 105762fef451SLois Curfman McInnes 10589b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 10599b94acceSBarry Smith 1060*5615d1e5SSatish Balay #undef __FUNC__ 1061*5615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 10629b94acceSBarry Smith /*@ 10639b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1064272ac6f2SLois Curfman McInnes of a nonlinear solver. 10659b94acceSBarry Smith 10669b94acceSBarry Smith Input Parameter: 10679b94acceSBarry Smith . snes - the SNES context 10688ddd3da0SLois Curfman McInnes . x - the solution vector 10699b94acceSBarry Smith 1070272ac6f2SLois Curfman McInnes Notes: 1071272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1072272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1073272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1074272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1075272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1076272ac6f2SLois Curfman McInnes 10779b94acceSBarry Smith .keywords: SNES, nonlinear, setup 10789b94acceSBarry Smith 10799b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 10809b94acceSBarry Smith @*/ 10818ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 10829b94acceSBarry Smith { 1083272ac6f2SLois Curfman McInnes int ierr, flg; 108477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 108577c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 10868ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 10879b94acceSBarry Smith 1088c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1089c1f60f51SBarry Smith /* 1090c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1091dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1092c1f60f51SBarry Smith */ 1093c1f60f51SBarry Smith if (flg) { 1094c1f60f51SBarry Smith Mat J; 1095c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1096c1f60f51SBarry Smith PLogObjectParent(snes,J); 1097c1f60f51SBarry Smith snes->mfshell = J; 1098c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1099c1f60f51SBarry Smith snes->jacobian = J; 1100c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1101c1f60f51SBarry Smith } 1102c1f60f51SBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1103c1f60f51SBarry Smith snes->jacobian = J; 1104c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1105e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free operator option"); 1106c1f60f51SBarry Smith } 1107272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1108c1f60f51SBarry Smith /* 1109dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1110c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1111c1f60f51SBarry Smith */ 111231615d04SLois Curfman McInnes if (flg) { 1113272ac6f2SLois Curfman McInnes Mat J; 1114272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1115272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1116272ac6f2SLois Curfman McInnes snes->mfshell = J; 111783e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 111883e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 111994a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 112083e56afdSLois Curfman McInnes } 112183e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 112283e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 112394a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1124e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free option"); 1125272ac6f2SLois Curfman McInnes } 11269b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1127e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1128e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1129e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first"); 1130e3372554SBarry Smith if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector"); 1131a703fe33SLois Curfman McInnes 1132a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 113340191667SLois Curfman McInnes if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1134a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1135a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1136a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1137a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1138a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1139a703fe33SLois Curfman McInnes } 11409b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1141e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first"); 1142e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first"); 114337fcc0dbSBarry Smith if (!snes->computeumfunction) 1144e3372554SBarry Smith SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first"); 1145e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first"); 1146e3372554SBarry Smith } else SETERRQ(1,0,"Unknown method class"); 1147a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1148a703fe33SLois Curfman McInnes snes->setup_called = 1; 1149a703fe33SLois Curfman McInnes return 0; 11509b94acceSBarry Smith } 11519b94acceSBarry Smith 1152*5615d1e5SSatish Balay #undef __FUNC__ 1153*5615d1e5SSatish Balay #define __FUNC__ "SNESDestroy" 11549b94acceSBarry Smith /*@C 11559b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 11569b94acceSBarry Smith with SNESCreate(). 11579b94acceSBarry Smith 11589b94acceSBarry Smith Input Parameter: 11599b94acceSBarry Smith . snes - the SNES context 11609b94acceSBarry Smith 11619b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 11629b94acceSBarry Smith 116363a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 11649b94acceSBarry Smith @*/ 11659b94acceSBarry Smith int SNESDestroy(SNES snes) 11669b94acceSBarry Smith { 11679b94acceSBarry Smith int ierr; 116877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 11699b94acceSBarry Smith ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr); 11700452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 11719b94acceSBarry Smith if (snes->mfshell) MatDestroy(snes->mfshell); 11729b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 11733e2e452bSBarry Smith if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 11746f24a144SLois Curfman McInnes if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 11759b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 11760452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 11779b94acceSBarry Smith return 0; 11789b94acceSBarry Smith } 11799b94acceSBarry Smith 11809b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 11819b94acceSBarry Smith 1182*5615d1e5SSatish Balay #undef __FUNC__ 1183*5615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 11849b94acceSBarry Smith /*@ 1185d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 11869b94acceSBarry Smith 11879b94acceSBarry Smith Input Parameters: 11889b94acceSBarry Smith . snes - the SNES context 118933174efeSLois Curfman McInnes . atol - absolute convergence tolerance 119033174efeSLois Curfman McInnes . rtol - relative convergence tolerance 119133174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 119233174efeSLois Curfman McInnes of the change in the solution between steps 119333174efeSLois Curfman McInnes . maxit - maximum number of iterations 119433174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 11959b94acceSBarry Smith 119633174efeSLois Curfman McInnes Options Database Keys: 119733174efeSLois Curfman McInnes $ -snes_atol <atol> 119833174efeSLois Curfman McInnes $ -snes_rtol <rtol> 119933174efeSLois Curfman McInnes $ -snes_stol <stol> 120033174efeSLois Curfman McInnes $ -snes_max_it <maxit> 120133174efeSLois Curfman McInnes $ -snes_max_funcs <maxf> 12029b94acceSBarry Smith 1203d7a720efSLois Curfman McInnes Notes: 12049b94acceSBarry Smith The default maximum number of iterations is 50. 12059b94acceSBarry Smith The default maximum number of function evaluations is 1000. 12069b94acceSBarry Smith 120733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 12089b94acceSBarry Smith 1209d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 12109b94acceSBarry Smith @*/ 1211d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 12129b94acceSBarry Smith { 121377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1214d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1215d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1216d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1217d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1218d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 12199b94acceSBarry Smith return 0; 12209b94acceSBarry Smith } 12219b94acceSBarry Smith 1222*5615d1e5SSatish Balay #undef __FUNC__ 1223*5615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 12249b94acceSBarry Smith /*@ 122533174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 122633174efeSLois Curfman McInnes 122733174efeSLois Curfman McInnes Input Parameters: 122833174efeSLois Curfman McInnes . snes - the SNES context 122933174efeSLois Curfman McInnes . atol - absolute convergence tolerance 123033174efeSLois Curfman McInnes . rtol - relative convergence tolerance 123133174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 123233174efeSLois Curfman McInnes of the change in the solution between steps 123333174efeSLois Curfman McInnes . maxit - maximum number of iterations 123433174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 123533174efeSLois Curfman McInnes 123633174efeSLois Curfman McInnes Notes: 123733174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 123833174efeSLois Curfman McInnes 123933174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 124033174efeSLois Curfman McInnes 124133174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 124233174efeSLois Curfman McInnes @*/ 124333174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 124433174efeSLois Curfman McInnes { 124533174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 124633174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 124733174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 124833174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 124933174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 125033174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 125133174efeSLois Curfman McInnes return 0; 125233174efeSLois Curfman McInnes } 125333174efeSLois Curfman McInnes 1254*5615d1e5SSatish Balay #undef __FUNC__ 1255*5615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 125633174efeSLois Curfman McInnes /*@ 12579b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 12589b94acceSBarry Smith 12599b94acceSBarry Smith Input Parameters: 12609b94acceSBarry Smith . snes - the SNES context 12619b94acceSBarry Smith . tol - tolerance 12629b94acceSBarry Smith 12639b94acceSBarry Smith Options Database Key: 1264d7a720efSLois Curfman McInnes $ -snes_trtol <tol> 12659b94acceSBarry Smith 12669b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 12679b94acceSBarry Smith 1268d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 12699b94acceSBarry Smith @*/ 12709b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 12719b94acceSBarry Smith { 127277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12739b94acceSBarry Smith snes->deltatol = tol; 12749b94acceSBarry Smith return 0; 12759b94acceSBarry Smith } 12769b94acceSBarry Smith 1277*5615d1e5SSatish Balay #undef __FUNC__ 1278*5615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 12799b94acceSBarry Smith /*@ 128077c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 12819b94acceSBarry Smith for unconstrained minimization solvers. 12829b94acceSBarry Smith 12839b94acceSBarry Smith Input Parameters: 12849b94acceSBarry Smith . snes - the SNES context 12859b94acceSBarry Smith . ftol - minimum function tolerance 12869b94acceSBarry Smith 12879b94acceSBarry Smith Options Database Key: 1288d7a720efSLois Curfman McInnes $ -snes_fmin <ftol> 12899b94acceSBarry Smith 12909b94acceSBarry Smith Note: 129177c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 12929b94acceSBarry Smith methods only. 12939b94acceSBarry Smith 12949b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 12959b94acceSBarry Smith 1296d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 12979b94acceSBarry Smith @*/ 129877c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 12999b94acceSBarry Smith { 130077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13019b94acceSBarry Smith snes->fmin = ftol; 13029b94acceSBarry Smith return 0; 13039b94acceSBarry Smith } 13049b94acceSBarry Smith 13059b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 13069b94acceSBarry Smith 1307*5615d1e5SSatish Balay #undef __FUNC__ 1308*5615d1e5SSatish Balay #define __FUNC__ "SNESSetMonitor" 13099b94acceSBarry Smith /*@C 13109b94acceSBarry Smith SNESSetMonitor - Sets the function that is to be used at every 13119b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 13129b94acceSBarry Smith progress. 13139b94acceSBarry Smith 13149b94acceSBarry Smith Input Parameters: 13159b94acceSBarry Smith . snes - the SNES context 13169b94acceSBarry Smith . func - monitoring routine 1317044dda88SLois Curfman McInnes . mctx - [optional] user-defined context for private data for the 1318044dda88SLois Curfman McInnes monitor routine (may be PETSC_NULL) 13199b94acceSBarry Smith 13209b94acceSBarry Smith Calling sequence of func: 1321682d7d0cSBarry Smith int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 13229b94acceSBarry Smith 13239b94acceSBarry Smith $ snes - the SNES context 13249b94acceSBarry Smith $ its - iteration number 1325044dda88SLois Curfman McInnes $ mctx - [optional] monitoring context 13269b94acceSBarry Smith $ 13279b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13289b94acceSBarry Smith $ norm - 2-norm function value (may be estimated) 13299b94acceSBarry Smith $ 13309b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13319b94acceSBarry Smith $ norm - 2-norm gradient value (may be estimated) 13329b94acceSBarry Smith 13339665c990SLois Curfman McInnes Options Database Keys: 13349665c990SLois Curfman McInnes $ -snes_monitor : sets SNESDefaultMonitor() 13359665c990SLois Curfman McInnes $ -snes_xmonitor : sets line graph monitor, 13369665c990SLois Curfman McInnes $ uses SNESLGMonitorCreate() 13379665c990SLois Curfman McInnes $ -snes_cancelmonitors : cancels all monitors that have 13389665c990SLois Curfman McInnes $ been hardwired into a code by 13399665c990SLois Curfman McInnes $ calls to SNESSetMonitor(), but 13409665c990SLois Curfman McInnes $ does not cancel those set via 13419665c990SLois Curfman McInnes $ the options database. 13429665c990SLois Curfman McInnes 13439665c990SLois Curfman McInnes 1344639f9d9dSBarry Smith Notes: 13456bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 13466bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 13476bc08f3fSLois Curfman McInnes order in which they were set. 1348639f9d9dSBarry Smith 13499b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 13509b94acceSBarry Smith 13519b94acceSBarry Smith .seealso: SNESDefaultMonitor() 13529b94acceSBarry Smith @*/ 135374679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 13549b94acceSBarry Smith { 1355639f9d9dSBarry Smith if (!func) { 1356639f9d9dSBarry Smith snes->numbermonitors = 0; 1357639f9d9dSBarry Smith return 0; 1358639f9d9dSBarry Smith } 1359639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1360e3372554SBarry Smith SETERRQ(1,0,"Too many monitors set"); 1361639f9d9dSBarry Smith } 1362639f9d9dSBarry Smith 1363639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1364639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 13659b94acceSBarry Smith return 0; 13669b94acceSBarry Smith } 13679b94acceSBarry Smith 1368*5615d1e5SSatish Balay #undef __FUNC__ 1369*5615d1e5SSatish Balay #define __FUNC__ "SNESSetConvergenceTest" 13709b94acceSBarry Smith /*@C 13719b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 13729b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 13739b94acceSBarry Smith 13749b94acceSBarry Smith Input Parameters: 13759b94acceSBarry Smith . snes - the SNES context 13769b94acceSBarry Smith . func - routine to test for convergence 1377044dda88SLois Curfman McInnes . cctx - [optional] context for private data for the convergence routine 1378044dda88SLois Curfman McInnes (may be PETSC_NULL) 13799b94acceSBarry Smith 13809b94acceSBarry Smith Calling sequence of func: 13819b94acceSBarry Smith int func (SNES snes,double xnorm,double gnorm, 13829b94acceSBarry Smith double f,void *cctx) 13839b94acceSBarry Smith 13849b94acceSBarry Smith $ snes - the SNES context 1385044dda88SLois Curfman McInnes $ cctx - [optional] convergence context 13869b94acceSBarry Smith $ xnorm - 2-norm of current iterate 13879b94acceSBarry Smith $ 13889b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13899b94acceSBarry Smith $ gnorm - 2-norm of current step 13909b94acceSBarry Smith $ f - 2-norm of function 13919b94acceSBarry Smith $ 13929b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13939b94acceSBarry Smith $ gnorm - 2-norm of current gradient 13949b94acceSBarry Smith $ f - function value 13959b94acceSBarry Smith 13969b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 13979b94acceSBarry Smith 139840191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 139940191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 14009b94acceSBarry Smith @*/ 140174679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 14029b94acceSBarry Smith { 14039b94acceSBarry Smith (snes)->converged = func; 14049b94acceSBarry Smith (snes)->cnvP = cctx; 14059b94acceSBarry Smith return 0; 14069b94acceSBarry Smith } 14079b94acceSBarry Smith 1408*5615d1e5SSatish Balay #undef __FUNC__ 1409*5615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 14109b94acceSBarry Smith /* 14119b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 14129b94acceSBarry Smith positive parameter delta. 14139b94acceSBarry Smith 14149b94acceSBarry Smith Input Parameters: 14159b94acceSBarry Smith . snes - the SNES context 14169b94acceSBarry Smith . y - approximate solution of linear system 14179b94acceSBarry Smith . fnorm - 2-norm of current function 14189b94acceSBarry Smith . delta - trust region size 14199b94acceSBarry Smith 14209b94acceSBarry Smith Output Parameters: 14219b94acceSBarry Smith . gpnorm - predicted function norm at the new point, assuming local 14229b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 14239b94acceSBarry Smith region, and exceeds zero otherwise. 14249b94acceSBarry Smith . ynorm - 2-norm of the step 14259b94acceSBarry Smith 14269b94acceSBarry Smith Note: 142740191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 14289b94acceSBarry Smith is set to be the maximum allowable step size. 14299b94acceSBarry Smith 14309b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 14319b94acceSBarry Smith */ 14329b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 14339b94acceSBarry Smith double *gpnorm,double *ynorm) 14349b94acceSBarry Smith { 14359b94acceSBarry Smith double norm; 14369b94acceSBarry Smith Scalar cnorm; 1437cddf8d76SBarry Smith VecNorm(y,NORM_2, &norm ); 14389b94acceSBarry Smith if (norm > *delta) { 14399b94acceSBarry Smith norm = *delta/norm; 14409b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 14419b94acceSBarry Smith cnorm = norm; 14429b94acceSBarry Smith VecScale( &cnorm, y ); 14439b94acceSBarry Smith *ynorm = *delta; 14449b94acceSBarry Smith } else { 14459b94acceSBarry Smith *gpnorm = 0.0; 14469b94acceSBarry Smith *ynorm = norm; 14479b94acceSBarry Smith } 14489b94acceSBarry Smith return 0; 14499b94acceSBarry Smith } 14509b94acceSBarry Smith 1451*5615d1e5SSatish Balay #undef __FUNC__ 1452*5615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 14539b94acceSBarry Smith /*@ 14549b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 145563a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 14569b94acceSBarry Smith 14579b94acceSBarry Smith Input Parameter: 14589b94acceSBarry Smith . snes - the SNES context 14598ddd3da0SLois Curfman McInnes . x - the solution vector 14609b94acceSBarry Smith 14619b94acceSBarry Smith Output Parameter: 14629b94acceSBarry Smith its - number of iterations until termination 14639b94acceSBarry Smith 14648ddd3da0SLois Curfman McInnes Note: 14658ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 14668ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 14678ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 14688ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 14698ddd3da0SLois Curfman McInnes 14709b94acceSBarry Smith .keywords: SNES, nonlinear, solve 14719b94acceSBarry Smith 147263a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 14739b94acceSBarry Smith @*/ 14748ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 14759b94acceSBarry Smith { 14763c7409f5SSatish Balay int ierr, flg; 1477052efed2SBarry Smith 147877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 147974679c65SBarry Smith PetscValidIntPointer(its); 1480c4fc05e7SBarry Smith if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1481c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 14829b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 14839b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 14849b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 14853c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 14866d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 14879b94acceSBarry Smith return 0; 14889b94acceSBarry Smith } 14899b94acceSBarry Smith 14909b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 149137fcc0dbSBarry Smith static NRList *__SNESList = 0; 14929b94acceSBarry Smith 1493*5615d1e5SSatish Balay #undef __FUNC__ 1494*5615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 14959b94acceSBarry Smith /*@ 14964b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 14979b94acceSBarry Smith 14989b94acceSBarry Smith Input Parameters: 14999b94acceSBarry Smith . snes - the SNES context 15009b94acceSBarry Smith . method - a known method 15019b94acceSBarry Smith 1502ae12b187SLois Curfman McInnes Options Database Command: 1503ae12b187SLois Curfman McInnes $ -snes_type <method> 1504ae12b187SLois Curfman McInnes $ Use -help for a list of available methods 1505ae12b187SLois Curfman McInnes $ (for instance, ls or tr) 1506ae12b187SLois Curfman McInnes 15079b94acceSBarry Smith Notes: 15089b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 15099b94acceSBarry Smith $ Systems of nonlinear equations: 151040191667SLois Curfman McInnes $ SNES_EQ_LS - Newton's method with line search 151140191667SLois Curfman McInnes $ SNES_EQ_TR - Newton's method with trust region 15129b94acceSBarry Smith $ Unconstrained minimization: 151340191667SLois Curfman McInnes $ SNES_UM_TR - Newton's method with trust region 151440191667SLois Curfman McInnes $ SNES_UM_LS - Newton's method with line search 15159b94acceSBarry Smith 1516ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1517ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1518ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1519ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1520ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1521ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1522ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1523ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1524ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1525ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1526a703fe33SLois Curfman McInnes 1527f536c699SSatish Balay .keywords: SNES, set, method 15289b94acceSBarry Smith @*/ 15294b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 15309b94acceSBarry Smith { 15319b94acceSBarry Smith int (*r)(SNES); 153277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 15337d1a2b2bSBarry Smith if (snes->type == (int) method) return 0; 15347d1a2b2bSBarry Smith 15359b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 153637fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 1537e3372554SBarry Smith if (!__SNESList) {SETERRQ(1,0,"Could not get methods");} 153837fcc0dbSBarry Smith r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1539e3372554SBarry Smith if (!r) {SETERRQ(1,0,"Unknown method");} 15400452661fSBarry Smith if (snes->data) PetscFree(snes->data); 15419b94acceSBarry Smith snes->set_method_called = 1; 15429b94acceSBarry Smith return (*r)(snes); 15439b94acceSBarry Smith } 15449b94acceSBarry Smith 15459b94acceSBarry Smith /* --------------------------------------------------------------------- */ 1546*5615d1e5SSatish Balay #undef __FUNC__ 1547*5615d1e5SSatish Balay #define __FUNC__ "SNESRegister" 15489b94acceSBarry Smith /*@C 15499b94acceSBarry Smith SNESRegister - Adds the method to the nonlinear solver package, given 15504b0e389bSBarry Smith a function pointer and a nonlinear solver name of the type SNESType. 15519b94acceSBarry Smith 15529b94acceSBarry Smith Input Parameters: 155340191667SLois Curfman McInnes . name - for instance SNES_EQ_LS, SNES_EQ_TR, ... 155440191667SLois Curfman McInnes . sname - corresponding string for name 15559b94acceSBarry Smith . create - routine to create method context 15569b94acceSBarry Smith 15579b94acceSBarry Smith .keywords: SNES, nonlinear, register 15589b94acceSBarry Smith 15599b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy() 15609b94acceSBarry Smith @*/ 15619b94acceSBarry Smith int SNESRegister(int name, char *sname, int (*create)(SNES)) 15629b94acceSBarry Smith { 15639b94acceSBarry Smith int ierr; 156437fcc0dbSBarry Smith if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 156537fcc0dbSBarry Smith NRRegister( __SNESList, name, sname, (int (*)(void*))create ); 15669b94acceSBarry Smith return 0; 15679b94acceSBarry Smith } 1568a847f771SSatish Balay 15699b94acceSBarry Smith /* --------------------------------------------------------------------- */ 1570*5615d1e5SSatish Balay #undef __FUNC__ 1571*5615d1e5SSatish Balay #define __FUNC__ "SNESRegisterDestroy" 15729b94acceSBarry Smith /*@C 15739b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 15749b94acceSBarry Smith registered by SNESRegister(). 15759b94acceSBarry Smith 15769b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 15779b94acceSBarry Smith 15789b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 15799b94acceSBarry Smith @*/ 15809b94acceSBarry Smith int SNESRegisterDestroy() 15819b94acceSBarry Smith { 158237fcc0dbSBarry Smith if (__SNESList) { 158337fcc0dbSBarry Smith NRDestroy( __SNESList ); 158437fcc0dbSBarry Smith __SNESList = 0; 15859b94acceSBarry Smith } 15869b94acceSBarry Smith return 0; 15879b94acceSBarry Smith } 15889b94acceSBarry Smith 1589*5615d1e5SSatish Balay #undef __FUNC__ 1590*5615d1e5SSatish Balay #define __FUNC__ "SNESGetTypeFromOptions_Private" 15919b94acceSBarry Smith /* 15924b0e389bSBarry Smith SNESGetTypeFromOptions_Private - Sets the selected method from the 15939b94acceSBarry Smith options database. 15949b94acceSBarry Smith 15959b94acceSBarry Smith Input Parameter: 15969b94acceSBarry Smith . ctx - the SNES context 15979b94acceSBarry Smith 15989b94acceSBarry Smith Output Parameter: 15999b94acceSBarry Smith . method - solver method 16009b94acceSBarry Smith 16019b94acceSBarry Smith Returns: 16029b94acceSBarry Smith Returns 1 if the method is found; 0 otherwise. 16039b94acceSBarry Smith 16049b94acceSBarry Smith Options Database Key: 16054b0e389bSBarry Smith $ -snes_type method 16069b94acceSBarry Smith */ 1607052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 16089b94acceSBarry Smith { 1609052efed2SBarry Smith int ierr; 16109b94acceSBarry Smith char sbuf[50]; 16115baf8537SBarry Smith 1612052efed2SBarry Smith ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1613052efed2SBarry Smith if (*flg) { 1614052efed2SBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 161537fcc0dbSBarry Smith *method = (SNESType)NRFindID( __SNESList, sbuf ); 16169b94acceSBarry Smith } 16179b94acceSBarry Smith return 0; 16189b94acceSBarry Smith } 16199b94acceSBarry Smith 1620*5615d1e5SSatish Balay #undef __FUNC__ 1621*5615d1e5SSatish Balay #define __FUNC__ "SNESGetType" 16229b94acceSBarry Smith /*@C 16239a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 16249b94acceSBarry Smith 16259b94acceSBarry Smith Input Parameter: 16264b0e389bSBarry Smith . snes - nonlinear solver context 16279b94acceSBarry Smith 16289b94acceSBarry Smith Output Parameter: 16299a28b0a6SLois Curfman McInnes . method - SNES method (or use PETSC_NULL) 16309a28b0a6SLois Curfman McInnes . name - name of SNES method (or use PETSC_NULL) 16319b94acceSBarry Smith 16329b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 16339b94acceSBarry Smith @*/ 16344b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name) 16359b94acceSBarry Smith { 163606a9b0adSLois Curfman McInnes int ierr; 163737fcc0dbSBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 16384b0e389bSBarry Smith if (method) *method = (SNESType) snes->type; 163937fcc0dbSBarry Smith if (name) *name = NRFindName( __SNESList, (int) snes->type ); 16409b94acceSBarry Smith return 0; 16419b94acceSBarry Smith } 16429b94acceSBarry Smith 16439b94acceSBarry Smith #include <stdio.h> 1644*5615d1e5SSatish Balay #undef __FUNC__ 1645*5615d1e5SSatish Balay #define __FUNC__ "SNESPrintTypes_Private" 16469b94acceSBarry Smith /* 16474b0e389bSBarry Smith SNESPrintTypes_Private - Prints the SNES methods available from the 16489b94acceSBarry Smith options database. 16499b94acceSBarry Smith 16509b94acceSBarry Smith Input Parameters: 165133455573SSatish Balay . comm - communicator (usually MPI_COMM_WORLD) 16529b94acceSBarry Smith . prefix - prefix (usually "-") 16534b0e389bSBarry Smith . name - the options database name (by default "snes_type") 16549b94acceSBarry Smith */ 165533455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 16569b94acceSBarry Smith { 16579b94acceSBarry Smith FuncList *entry; 165837fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 165937fcc0dbSBarry Smith entry = __SNESList->head; 166077c4ece6SBarry Smith PetscPrintf(comm," %s%s (one of)",prefix,name); 16619b94acceSBarry Smith while (entry) { 166277c4ece6SBarry Smith PetscPrintf(comm," %s",entry->name); 16639b94acceSBarry Smith entry = entry->next; 16649b94acceSBarry Smith } 166577c4ece6SBarry Smith PetscPrintf(comm,"\n"); 16669b94acceSBarry Smith return 0; 16679b94acceSBarry Smith } 16689b94acceSBarry Smith 1669*5615d1e5SSatish Balay #undef __FUNC__ 1670*5615d1e5SSatish Balay #define __FUNC__ "SNESGetSolution" 16719b94acceSBarry Smith /*@C 16729b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 16739b94acceSBarry Smith stored. 16749b94acceSBarry Smith 16759b94acceSBarry Smith Input Parameter: 16769b94acceSBarry Smith . snes - the SNES context 16779b94acceSBarry Smith 16789b94acceSBarry Smith Output Parameter: 16799b94acceSBarry Smith . x - the solution 16809b94acceSBarry Smith 16819b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 16829b94acceSBarry Smith 1683a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 16849b94acceSBarry Smith @*/ 16859b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 16869b94acceSBarry Smith { 168777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16889b94acceSBarry Smith *x = snes->vec_sol_always; 16899b94acceSBarry Smith return 0; 16909b94acceSBarry Smith } 16919b94acceSBarry Smith 1692*5615d1e5SSatish Balay #undef __FUNC__ 1693*5615d1e5SSatish Balay #define __FUNC__ "SNESGetSolutionUpdate" 16949b94acceSBarry Smith /*@C 16959b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 16969b94acceSBarry Smith stored. 16979b94acceSBarry Smith 16989b94acceSBarry Smith Input Parameter: 16999b94acceSBarry Smith . snes - the SNES context 17009b94acceSBarry Smith 17019b94acceSBarry Smith Output Parameter: 17029b94acceSBarry Smith . x - the solution update 17039b94acceSBarry Smith 17049b94acceSBarry Smith Notes: 17059b94acceSBarry Smith This vector is implementation dependent. 17069b94acceSBarry Smith 17079b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 17089b94acceSBarry Smith 17099b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 17109b94acceSBarry Smith @*/ 17119b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 17129b94acceSBarry Smith { 171377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17149b94acceSBarry Smith *x = snes->vec_sol_update_always; 17159b94acceSBarry Smith return 0; 17169b94acceSBarry Smith } 17179b94acceSBarry Smith 1718*5615d1e5SSatish Balay #undef __FUNC__ 1719*5615d1e5SSatish Balay #define __FUNC__ "SNESGetFunction" 17209b94acceSBarry Smith /*@C 17213638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 17229b94acceSBarry Smith 17239b94acceSBarry Smith Input Parameter: 17249b94acceSBarry Smith . snes - the SNES context 17259b94acceSBarry Smith 17269b94acceSBarry Smith Output Parameter: 17273638b69dSLois Curfman McInnes . r - the function 17289b94acceSBarry Smith 17299b94acceSBarry Smith Notes: 17309b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 17319b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 17329b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 17339b94acceSBarry Smith 1734a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 17359b94acceSBarry Smith 173631615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 173731615d04SLois Curfman McInnes SNESGetGradient() 17389b94acceSBarry Smith @*/ 17399b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 17409b94acceSBarry Smith { 174177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1742e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0, 17437c792091SSatish Balay "For SNES_NONLINEAR_EQUATIONS only"); 17449b94acceSBarry Smith *r = snes->vec_func_always; 17459b94acceSBarry Smith return 0; 17469b94acceSBarry Smith } 17479b94acceSBarry Smith 1748*5615d1e5SSatish Balay #undef __FUNC__ 1749*5615d1e5SSatish Balay #define __FUNC__ "SNESGetGradient" 17509b94acceSBarry Smith /*@C 17513638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 17529b94acceSBarry Smith 17539b94acceSBarry Smith Input Parameter: 17549b94acceSBarry Smith . snes - the SNES context 17559b94acceSBarry Smith 17569b94acceSBarry Smith Output Parameter: 17579b94acceSBarry Smith . r - the gradient 17589b94acceSBarry Smith 17599b94acceSBarry Smith Notes: 17609b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 17619b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 17629b94acceSBarry Smith SNESGetFunction(). 17639b94acceSBarry Smith 17649b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 17659b94acceSBarry Smith 176631615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 17679b94acceSBarry Smith @*/ 17689b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 17699b94acceSBarry Smith { 177077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1771e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 177263c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 17739b94acceSBarry Smith *r = snes->vec_func_always; 17749b94acceSBarry Smith return 0; 17759b94acceSBarry Smith } 17769b94acceSBarry Smith 1777*5615d1e5SSatish Balay #undef __FUNC__ 1778*5615d1e5SSatish Balay #define __FUNC__ "SNESGetMinimizationFunction" 17799b94acceSBarry Smith /*@ 17809b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 17819b94acceSBarry Smith unconstrained minimization problems. 17829b94acceSBarry Smith 17839b94acceSBarry Smith Input Parameter: 17849b94acceSBarry Smith . snes - the SNES context 17859b94acceSBarry Smith 17869b94acceSBarry Smith Output Parameter: 17879b94acceSBarry Smith . r - the function 17889b94acceSBarry Smith 17899b94acceSBarry Smith Notes: 17909b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 17919b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 17929b94acceSBarry Smith SNESGetFunction(). 17939b94acceSBarry Smith 17949b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 17959b94acceSBarry Smith 179631615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 17979b94acceSBarry Smith @*/ 17989b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 17999b94acceSBarry Smith { 180077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 180174679c65SBarry Smith PetscValidScalarPointer(r); 1802e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 180363c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 18049b94acceSBarry Smith *r = snes->fc; 18059b94acceSBarry Smith return 0; 18069b94acceSBarry Smith } 18079b94acceSBarry Smith 1808*5615d1e5SSatish Balay #undef __FUNC__ 1809*5615d1e5SSatish Balay #define __FUNC__ "SNESSetOptionsPrefix" 18103c7409f5SSatish Balay /*@C 18113c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1812639f9d9dSBarry Smith SNES options in the database. You must NOT include the - at the beginning of 1813639f9d9dSBarry Smith the prefix name. 18143c7409f5SSatish Balay 18153c7409f5SSatish Balay Input Parameter: 18163c7409f5SSatish Balay . snes - the SNES context 18173c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 18183c7409f5SSatish Balay 18193c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1820a86d99e1SLois Curfman McInnes 1821a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 18223c7409f5SSatish Balay @*/ 18233c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 18243c7409f5SSatish Balay { 18253c7409f5SSatish Balay int ierr; 18263c7409f5SSatish Balay 182777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1828639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 18293c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 18303c7409f5SSatish Balay return 0; 18313c7409f5SSatish Balay } 18323c7409f5SSatish Balay 1833*5615d1e5SSatish Balay #undef __FUNC__ 1834*5615d1e5SSatish Balay #define __FUNC__ "SNESAppendOptionsPrefix" 18353c7409f5SSatish Balay /*@C 1836f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1837639f9d9dSBarry Smith SNES options in the database. You must NOT include the - at the beginning of 1838639f9d9dSBarry Smith the prefix name. 18393c7409f5SSatish Balay 18403c7409f5SSatish Balay Input Parameter: 18413c7409f5SSatish Balay . snes - the SNES context 18423c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 18433c7409f5SSatish Balay 18443c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1845a86d99e1SLois Curfman McInnes 1846a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 18473c7409f5SSatish Balay @*/ 18483c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 18493c7409f5SSatish Balay { 18503c7409f5SSatish Balay int ierr; 18513c7409f5SSatish Balay 185277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1853639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 18543c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 18553c7409f5SSatish Balay return 0; 18563c7409f5SSatish Balay } 18573c7409f5SSatish Balay 1858*5615d1e5SSatish Balay #undef __FUNC__ 1859*5615d1e5SSatish Balay #define __FUNC__ "SNESGetOptionsPrefix" 1860c83590e2SSatish Balay /*@ 18613c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 18623c7409f5SSatish Balay SNES options in the database. 18633c7409f5SSatish Balay 18643c7409f5SSatish Balay Input Parameter: 18653c7409f5SSatish Balay . snes - the SNES context 18663c7409f5SSatish Balay 18673c7409f5SSatish Balay Output Parameter: 18683c7409f5SSatish Balay . prefix - pointer to the prefix string used 18693c7409f5SSatish Balay 18703c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 1871a86d99e1SLois Curfman McInnes 1872a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 18733c7409f5SSatish Balay @*/ 18743c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 18753c7409f5SSatish Balay { 18763c7409f5SSatish Balay int ierr; 18773c7409f5SSatish Balay 187877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1879639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 18803c7409f5SSatish Balay return 0; 18813c7409f5SSatish Balay } 18823c7409f5SSatish Balay 18833c7409f5SSatish Balay 18849b94acceSBarry Smith 18859b94acceSBarry Smith 18869b94acceSBarry Smith 1887