19b94acceSBarry Smith #ifndef lint 2*93c39befSBarry Smith static char vcid[] = "$Id: snes.c,v 1.123 1997/06/05 12:57:14 bsmith Exp bsmith $"; 39b94acceSBarry Smith #endif 49b94acceSBarry Smith 570f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 6f5eb4b81SSatish Balay #include "src/sys/nreg.h" 79b94acceSBarry Smith #include "pinclude/pviewer.h" 89b94acceSBarry Smith #include <math.h> 99b94acceSBarry Smith 10052efed2SBarry Smith extern int SNESGetTypeFromOptions_Private(SNES,SNESType*,int*); 1133455573SSatish Balay extern int SNESPrintTypes_Private(MPI_Comm,char*,char*); 129b94acceSBarry Smith 1384cb2905SBarry Smith int SNESRegisterAllCalled = 0; 1484cb2905SBarry Smith 155615d1e5SSatish Balay #undef __FUNC__ 165eea60f9SBarry Smith #define __FUNC__ "SNESView" /* ADIC Ignore */ 179b94acceSBarry Smith /*@ 189b94acceSBarry Smith SNESView - Prints the SNES data structure. 199b94acceSBarry Smith 209b94acceSBarry Smith Input Parameters: 219b94acceSBarry Smith . SNES - the SNES context 229b94acceSBarry Smith . viewer - visualization context 239b94acceSBarry Smith 249b94acceSBarry Smith Options Database Key: 259b94acceSBarry Smith $ -snes_view : calls SNESView() at end of SNESSolve() 269b94acceSBarry Smith 279b94acceSBarry Smith Notes: 289b94acceSBarry Smith The available visualization contexts include 296d4a8577SBarry Smith $ VIEWER_STDOUT_SELF - standard output (default) 306d4a8577SBarry Smith $ VIEWER_STDOUT_WORLD - synchronized standard 319b94acceSBarry Smith $ output where only the first processor opens 329b94acceSBarry Smith $ the file. All other processors send their 339b94acceSBarry Smith $ data to the first processor to print. 349b94acceSBarry Smith 359b94acceSBarry Smith The user can open alternative vistualization contexts with 36dbb450caSBarry Smith $ ViewerFileOpenASCII() - output to a specified file 379b94acceSBarry Smith 389b94acceSBarry Smith .keywords: SNES, view 399b94acceSBarry Smith 40dbb450caSBarry Smith .seealso: ViewerFileOpenASCII() 419b94acceSBarry Smith @*/ 429b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 439b94acceSBarry Smith { 449b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 459b94acceSBarry Smith FILE *fd; 469b94acceSBarry Smith int ierr; 479b94acceSBarry Smith SLES sles; 489b94acceSBarry Smith char *method; 4919bcc07fSBarry Smith ViewerType vtype; 509b94acceSBarry Smith 5174679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5274679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5374679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5474679c65SBarry Smith 5519bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 5619bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 5790ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 5877c4ece6SBarry Smith PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 594b0e389bSBarry Smith SNESGetType(snes,PETSC_NULL,&method); 6077c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," method: %s\n",method); 619b94acceSBarry Smith if (snes->view) (*snes->view)((PetscObject)snes,viewer); 6277c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 639b94acceSBarry Smith " maximum iterations=%d, maximum function evaluations=%d\n", 649b94acceSBarry Smith snes->max_its,snes->max_funcs); 6577c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 669b94acceSBarry Smith " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 679b94acceSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol); 687a00f4a9SLois Curfman McInnes PetscFPrintf(snes->comm,fd, 697a00f4a9SLois Curfman McInnes " total number of linear solver iterations=%d\n",snes->linear_its); 70ae3c334cSLois Curfman McInnes PetscFPrintf(snes->comm,fd, 7168632166SLois Curfman McInnes " total number of function evaluations=%d\n",snes->nfuncs); 729b94acceSBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 7377c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 749b94acceSBarry Smith if (snes->ksp_ewconv) { 759b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 769b94acceSBarry Smith if (kctx) { 7777c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 789b94acceSBarry Smith " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 799b94acceSBarry Smith kctx->version); 8077c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 819b94acceSBarry Smith " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 829b94acceSBarry Smith kctx->rtol_max,kctx->threshold); 8377c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 849b94acceSBarry Smith kctx->gamma,kctx->alpha,kctx->alpha2); 859b94acceSBarry Smith } 869b94acceSBarry Smith } 87c30f285eSLois Curfman McInnes } else if (vtype == STRING_VIEWER) { 88c30f285eSLois Curfman McInnes SNESGetType(snes,PETSC_NULL,&method); 89c30f285eSLois Curfman McInnes ViewerStringSPrintf(viewer," %-3.3s",method); 9019bcc07fSBarry Smith } 919b94acceSBarry Smith SNESGetSLES(snes,&sles); 929b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 939b94acceSBarry Smith return 0; 949b94acceSBarry Smith } 959b94acceSBarry Smith 96639f9d9dSBarry Smith /* 97639f9d9dSBarry Smith We retain a list of functions that also take SNES command 98639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 99639f9d9dSBarry Smith */ 100639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 101639f9d9dSBarry Smith static int numberofsetfromoptions; 102639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 103639f9d9dSBarry Smith 1045615d1e5SSatish Balay #undef __FUNC__ 1055eea60f9SBarry Smith #define __FUNC__ "SNESAddOptionsChecker" /* ADIC Ignore */ 106639f9d9dSBarry Smith /*@ 107639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 108639f9d9dSBarry Smith 109639f9d9dSBarry Smith Input Parameter: 110639f9d9dSBarry Smith . snescheck - function that checks for options 111639f9d9dSBarry Smith 112639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 113639f9d9dSBarry Smith @*/ 114639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 115639f9d9dSBarry Smith { 116639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 117e3372554SBarry Smith SETERRQ(1,0,"Too many options checkers, only 5 allowed"); 118639f9d9dSBarry Smith } 119639f9d9dSBarry Smith 120639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 121639f9d9dSBarry Smith return 0; 122639f9d9dSBarry Smith } 123639f9d9dSBarry Smith 1245615d1e5SSatish Balay #undef __FUNC__ 1255615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1269b94acceSBarry Smith /*@ 127682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1289b94acceSBarry Smith 1299b94acceSBarry Smith Input Parameter: 1309b94acceSBarry Smith . snes - the SNES context 1319b94acceSBarry Smith 1329b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1339b94acceSBarry Smith 134a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1359b94acceSBarry Smith @*/ 1369b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1379b94acceSBarry Smith { 1384b0e389bSBarry Smith SNESType method; 1399b94acceSBarry Smith double tmp; 1409b94acceSBarry Smith SLES sles; 14181f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1423c7409f5SSatish Balay int version = PETSC_DEFAULT; 1439b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1449b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1459b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1469b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1479b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1489b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1499b94acceSBarry Smith 15081f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 15181f57ec7SBarry Smith 15281f57ec7SBarry Smith 15377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 154e3372554SBarry Smith if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp"); 155052efed2SBarry Smith ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 156052efed2SBarry Smith if (flg) { 157052efed2SBarry Smith ierr = SNESSetType(snes,method); CHKERRQ(ierr); 1589b94acceSBarry Smith } 1595334005bSBarry Smith else if (!snes->set_method_called) { 1605334005bSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 16140191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 1625334005bSBarry Smith } 1635334005bSBarry Smith else { 16440191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 1655334005bSBarry Smith } 1665334005bSBarry Smith } 1675334005bSBarry Smith 1683c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1693c7409f5SSatish Balay if (flg) { SNESPrintHelp(snes); } 1703c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 171d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);} 1723c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 173d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);} 1743c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 175d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); } 1763c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 1773c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 178d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 179d7a720efSLois Curfman McInnes if (flg) { SNESSetTrustRegionTolerance(snes,tmp); } 180d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 181d7a720efSLois Curfman McInnes if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); } 1823c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 1833c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 184b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 185b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 186b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 187b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 188b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 189b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 190b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 191*93c39befSBarry Smith 192*93c39befSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr); 193*93c39befSBarry Smith if (flg) {snes->converged = 0;} 194*93c39befSBarry Smith 1959b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 1969b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 1975bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 1985bbad29bSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);} 1993c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 200639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 2013c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 202d31a3109SLois Curfman McInnes if (flg) { 203639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 204d31a3109SLois Curfman McInnes } 20581f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2063c7409f5SSatish Balay if (flg) { 20717699dbbSLois Curfman McInnes int rank = 0; 208d7e8b826SBarry Smith DrawLG lg; 20917699dbbSLois Curfman McInnes MPI_Initialized(&rank); 21017699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 21117699dbbSLois Curfman McInnes if (!rank) { 21281f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2139b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 214f8c826e1SBarry Smith snes->xmonitor = lg; 2159b94acceSBarry Smith PLogObjectParent(snes,lg); 2169b94acceSBarry Smith } 2179b94acceSBarry Smith } 2183c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2193c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2209b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2219b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 22294a424c1SBarry Smith PLogInfo(snes, 22331615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 22431615d04SLois Curfman McInnes } 22531615d04SLois Curfman McInnes else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 22631615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 22731615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 22894a424c1SBarry Smith PLogInfo(snes, 22931615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2309b94acceSBarry Smith } 231639f9d9dSBarry Smith 232639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 233639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 234639f9d9dSBarry Smith } 235639f9d9dSBarry Smith 2369b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2379b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 2389b94acceSBarry Smith if (!snes->setfromoptions) return 0; 2399b94acceSBarry Smith return (*snes->setfromoptions)(snes); 2409b94acceSBarry Smith } 2419b94acceSBarry Smith 2425615d1e5SSatish Balay #undef __FUNC__ 2435eea60f9SBarry Smith #define __FUNC__ "SNESPrintHelp" /* ADIC Ignore */ 2449b94acceSBarry Smith /*@ 2459b94acceSBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2469b94acceSBarry Smith 2479b94acceSBarry Smith Input Parameter: 2489b94acceSBarry Smith . snes - the SNES context 2499b94acceSBarry Smith 250a703fe33SLois Curfman McInnes Options Database Keys: 251a703fe33SLois Curfman McInnes $ -help, -h 252a703fe33SLois Curfman McInnes 2539b94acceSBarry Smith .keywords: SNES, nonlinear, help 2549b94acceSBarry Smith 255682d7d0cSBarry Smith .seealso: SNESSetFromOptions() 2569b94acceSBarry Smith @*/ 2579b94acceSBarry Smith int SNESPrintHelp(SNES snes) 2589b94acceSBarry Smith { 2596daaf66cSBarry Smith char p[64]; 2606daaf66cSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2616daaf66cSBarry Smith 26277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2636daaf66cSBarry Smith 2646daaf66cSBarry Smith PetscStrcpy(p,"-"); 2656daaf66cSBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2666daaf66cSBarry Smith 2676daaf66cSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2686daaf66cSBarry Smith 269d31a3109SLois Curfman McInnes PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 27033455573SSatish Balay SNESPrintTypes_Private(snes->comm,p,"snes_type"); 27177c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 272b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 273b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 274b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 275b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 276b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 277b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 2785bbad29bSBarry Smith PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 2795bbad29bSBarry Smith PetscPrintf(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 280d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor\n",p); 281d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 282b18e04deSLois Curfman McInnes if (snes->type == SNES_NONLINEAR_EQUATIONS) { 28377c4ece6SBarry Smith PetscPrintf(snes->comm, 284d31a3109SLois Curfman McInnes " Options for solving systems of nonlinear equations only:\n"); 28577c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 28677c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 287ef1dfb11SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 288*93c39befSBarry Smith PetscPrintf(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p); 28977c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 29077c4ece6SBarry Smith PetscPrintf(snes->comm, 291b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 29277c4ece6SBarry Smith PetscPrintf(snes->comm, 293b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 29477c4ece6SBarry Smith PetscPrintf(snes->comm, 295b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 29677c4ece6SBarry Smith PetscPrintf(snes->comm, 297b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 29877c4ece6SBarry Smith PetscPrintf(snes->comm, 299b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 30077c4ece6SBarry Smith PetscPrintf(snes->comm, 301b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 30277c4ece6SBarry Smith PetscPrintf(snes->comm, 303b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 304b18e04deSLois Curfman McInnes } 305b18e04deSLois Curfman McInnes else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 30677c4ece6SBarry Smith PetscPrintf(snes->comm, 307d31a3109SLois Curfman McInnes " Options for solving unconstrained minimization problems only:\n"); 308b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 30977c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 31077c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 311b18e04deSLois Curfman McInnes } 3124537a946SLois Curfman McInnes PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 31377c4ece6SBarry Smith PetscPrintf(snes->comm,"a particular method\n"); 3146daaf66cSBarry Smith if (snes->printhelp) (*snes->printhelp)(snes,p); 3159b94acceSBarry Smith return 0; 3169b94acceSBarry Smith } 317a847f771SSatish Balay 3185615d1e5SSatish Balay #undef __FUNC__ 3195eea60f9SBarry Smith #define __FUNC__ "SNESSetApplicationContext" /* ADIC Ignore */ 3209b94acceSBarry Smith /*@ 3219b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 3229b94acceSBarry Smith the nonlinear solvers. 3239b94acceSBarry Smith 3249b94acceSBarry Smith Input Parameters: 3259b94acceSBarry Smith . snes - the SNES context 3269b94acceSBarry Smith . usrP - optional user context 3279b94acceSBarry Smith 3289b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3299b94acceSBarry Smith 3309b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3319b94acceSBarry Smith @*/ 3329b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3339b94acceSBarry Smith { 33477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3359b94acceSBarry Smith snes->user = usrP; 3369b94acceSBarry Smith return 0; 3379b94acceSBarry Smith } 33874679c65SBarry Smith 3395615d1e5SSatish Balay #undef __FUNC__ 3405eea60f9SBarry Smith #define __FUNC__ "SNESGetApplicationContext" /* ADIC Ignore */ 3419b94acceSBarry Smith /*@C 3429b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3439b94acceSBarry Smith nonlinear solvers. 3449b94acceSBarry Smith 3459b94acceSBarry Smith Input Parameter: 3469b94acceSBarry Smith . snes - SNES context 3479b94acceSBarry Smith 3489b94acceSBarry Smith Output Parameter: 3499b94acceSBarry Smith . usrP - user context 3509b94acceSBarry Smith 3519b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3529b94acceSBarry Smith 3539b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3549b94acceSBarry Smith @*/ 3559b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3569b94acceSBarry Smith { 35777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3589b94acceSBarry Smith *usrP = snes->user; 3599b94acceSBarry Smith return 0; 3609b94acceSBarry Smith } 36174679c65SBarry Smith 3625615d1e5SSatish Balay #undef __FUNC__ 3635eea60f9SBarry Smith #define __FUNC__ "SNESGetIterationNumber" /* ADIC Ignore */ 3649b94acceSBarry Smith /*@ 3659b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3669b94acceSBarry Smith nonlinear solver. 3679b94acceSBarry Smith 3689b94acceSBarry Smith Input Parameter: 3699b94acceSBarry Smith . snes - SNES context 3709b94acceSBarry Smith 3719b94acceSBarry Smith Output Parameter: 3729b94acceSBarry Smith . iter - iteration number 3739b94acceSBarry Smith 3749b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3759b94acceSBarry Smith @*/ 3769b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3779b94acceSBarry Smith { 37877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 37974679c65SBarry Smith PetscValidIntPointer(iter); 3809b94acceSBarry Smith *iter = snes->iter; 3819b94acceSBarry Smith return 0; 3829b94acceSBarry Smith } 38374679c65SBarry Smith 3845615d1e5SSatish Balay #undef __FUNC__ 3855615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 3869b94acceSBarry Smith /*@ 3879b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3889b94acceSBarry Smith with SNESSSetFunction(). 3899b94acceSBarry Smith 3909b94acceSBarry Smith Input Parameter: 3919b94acceSBarry Smith . snes - SNES context 3929b94acceSBarry Smith 3939b94acceSBarry Smith Output Parameter: 3949b94acceSBarry Smith . fnorm - 2-norm of function 3959b94acceSBarry Smith 3969b94acceSBarry Smith Note: 3979b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 398a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 399a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 4009b94acceSBarry Smith 4019b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 402a86d99e1SLois Curfman McInnes 403a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 4049b94acceSBarry Smith @*/ 4059b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 4069b94acceSBarry Smith { 40777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 40874679c65SBarry Smith PetscValidScalarPointer(fnorm); 40974679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 410d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 41174679c65SBarry Smith } 4129b94acceSBarry Smith *fnorm = snes->norm; 4139b94acceSBarry Smith return 0; 4149b94acceSBarry Smith } 41574679c65SBarry Smith 4165615d1e5SSatish Balay #undef __FUNC__ 4175615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 4189b94acceSBarry Smith /*@ 4199b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 4209b94acceSBarry Smith with SNESSSetGradient(). 4219b94acceSBarry Smith 4229b94acceSBarry Smith Input Parameter: 4239b94acceSBarry Smith . snes - SNES context 4249b94acceSBarry Smith 4259b94acceSBarry Smith Output Parameter: 4269b94acceSBarry Smith . fnorm - 2-norm of gradient 4279b94acceSBarry Smith 4289b94acceSBarry Smith Note: 4299b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 430a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 431a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4329b94acceSBarry Smith 4339b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 434a86d99e1SLois Curfman McInnes 435a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4369b94acceSBarry Smith @*/ 4379b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4389b94acceSBarry Smith { 43977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 44074679c65SBarry Smith PetscValidScalarPointer(gnorm); 44174679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 442d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 44374679c65SBarry Smith } 4449b94acceSBarry Smith *gnorm = snes->norm; 4459b94acceSBarry Smith return 0; 4469b94acceSBarry Smith } 44774679c65SBarry Smith 4485615d1e5SSatish Balay #undef __FUNC__ 4495eea60f9SBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" /* ADIC Ignore */ 4509b94acceSBarry Smith /*@ 4519b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4529b94acceSBarry Smith attempted by the nonlinear solver. 4539b94acceSBarry Smith 4549b94acceSBarry Smith Input Parameter: 4559b94acceSBarry Smith . snes - SNES context 4569b94acceSBarry Smith 4579b94acceSBarry Smith Output Parameter: 4589b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4599b94acceSBarry Smith 460c96a6f78SLois Curfman McInnes Notes: 461c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 462c96a6f78SLois Curfman McInnes 4639b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4649b94acceSBarry Smith @*/ 4659b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4669b94acceSBarry Smith { 46777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 46874679c65SBarry Smith PetscValidIntPointer(nfails); 4699b94acceSBarry Smith *nfails = snes->nfailures; 4709b94acceSBarry Smith return 0; 4719b94acceSBarry Smith } 472a847f771SSatish Balay 4735615d1e5SSatish Balay #undef __FUNC__ 4745eea60f9SBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" /* ADIC Ignore */ 475c96a6f78SLois Curfman McInnes /*@ 476c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 477c96a6f78SLois Curfman McInnes used by the nonlinear solver. 478c96a6f78SLois Curfman McInnes 479c96a6f78SLois Curfman McInnes Input Parameter: 480c96a6f78SLois Curfman McInnes . snes - SNES context 481c96a6f78SLois Curfman McInnes 482c96a6f78SLois Curfman McInnes Output Parameter: 483c96a6f78SLois Curfman McInnes . lits - number of linear iterations 484c96a6f78SLois Curfman McInnes 485c96a6f78SLois Curfman McInnes Notes: 486c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 487c96a6f78SLois Curfman McInnes 488c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 489c96a6f78SLois Curfman McInnes @*/ 49086bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 491c96a6f78SLois Curfman McInnes { 492c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 493c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 494c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 495c96a6f78SLois Curfman McInnes return 0; 496c96a6f78SLois Curfman McInnes } 497c96a6f78SLois Curfman McInnes 498c96a6f78SLois Curfman McInnes #undef __FUNC__ 4995eea60f9SBarry Smith #define __FUNC__ "SNESGetSLES" /* ADIC Ignore */ 5009b94acceSBarry Smith /*@C 5019b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 5029b94acceSBarry Smith 5039b94acceSBarry Smith Input Parameter: 5049b94acceSBarry Smith . snes - the SNES context 5059b94acceSBarry Smith 5069b94acceSBarry Smith Output Parameter: 5079b94acceSBarry Smith . sles - the SLES context 5089b94acceSBarry Smith 5099b94acceSBarry Smith Notes: 5109b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 5119b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5129b94acceSBarry Smith KSP and PC contexts as well. 5139b94acceSBarry Smith 5149b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 5159b94acceSBarry Smith 5169b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 5179b94acceSBarry Smith @*/ 5189b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 5199b94acceSBarry Smith { 52077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5219b94acceSBarry Smith *sles = snes->sles; 5229b94acceSBarry Smith return 0; 5239b94acceSBarry Smith } 5249b94acceSBarry Smith /* -----------------------------------------------------------*/ 5255615d1e5SSatish Balay #undef __FUNC__ 5265615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 5279b94acceSBarry Smith /*@C 5289b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5299b94acceSBarry Smith 5309b94acceSBarry Smith Input Parameter: 5319b94acceSBarry Smith . comm - MPI communicator 5329b94acceSBarry Smith . type - type of method, one of 5339b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS 5349b94acceSBarry Smith $ (for systems of nonlinear equations) 5359b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION 5369b94acceSBarry Smith $ (for unconstrained minimization) 5379b94acceSBarry Smith 5389b94acceSBarry Smith Output Parameter: 5399b94acceSBarry Smith . outsnes - the new SNES context 5409b94acceSBarry Smith 541c1f60f51SBarry Smith Options Database Key: 54219bd6747SLois Curfman McInnes $ -snes_mf - use default matrix-free Jacobian-vector products, 54319bd6747SLois Curfman McInnes $ and no preconditioning matrix 54419bd6747SLois Curfman McInnes $ -snes_mf_operator - use default matrix-free Jacobian-vector 54519bd6747SLois Curfman McInnes $ products, and a user-provided preconditioning matrix 54619bd6747SLois Curfman McInnes $ as set by SNESSetJacobian() 54719bd6747SLois Curfman McInnes $ -snes_fd - use (slow!) finite differences to compute Jacobian 548c1f60f51SBarry Smith 5499b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5509b94acceSBarry Smith 55163a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 5529b94acceSBarry Smith @*/ 5534b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 5549b94acceSBarry Smith { 5559b94acceSBarry Smith int ierr; 5569b94acceSBarry Smith SNES snes; 5579b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 55837fcc0dbSBarry Smith 5599b94acceSBarry Smith *outsnes = 0; 5602a0acf01SLois Curfman McInnes if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 561d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 562f09e8eb9SSatish Balay PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm); 5639b94acceSBarry Smith PLogObjectCreate(snes); 5649b94acceSBarry Smith snes->max_its = 50; 5659750a799SBarry Smith snes->max_funcs = 10000; 5669b94acceSBarry Smith snes->norm = 0.0; 5675a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 568b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 569b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5709b94acceSBarry Smith snes->atol = 1.e-10; 5715a2d0531SBarry Smith } 572b4874afaSBarry Smith else { 573b4874afaSBarry Smith snes->rtol = 1.e-8; 574b4874afaSBarry Smith snes->ttol = 0.0; 575b4874afaSBarry Smith snes->atol = 1.e-50; 576b4874afaSBarry Smith } 5779b94acceSBarry Smith snes->xtol = 1.e-8; 578b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5799b94acceSBarry Smith snes->nfuncs = 0; 5809b94acceSBarry Smith snes->nfailures = 0; 5817a00f4a9SLois Curfman McInnes snes->linear_its = 0; 582639f9d9dSBarry Smith snes->numbermonitors = 0; 5839b94acceSBarry Smith snes->data = 0; 5849b94acceSBarry Smith snes->view = 0; 5859b94acceSBarry Smith snes->computeumfunction = 0; 5869b94acceSBarry Smith snes->umfunP = 0; 5879b94acceSBarry Smith snes->fc = 0; 5889b94acceSBarry Smith snes->deltatol = 1.e-12; 5899b94acceSBarry Smith snes->fmin = -1.e30; 5909b94acceSBarry Smith snes->method_class = type; 5919b94acceSBarry Smith snes->set_method_called = 0; 592a703fe33SLois Curfman McInnes snes->setup_called = 0; 5939b94acceSBarry Smith snes->ksp_ewconv = 0; 5947d1a2b2bSBarry Smith snes->type = -1; 5956f24a144SLois Curfman McInnes snes->xmonitor = 0; 5966f24a144SLois Curfman McInnes snes->vwork = 0; 5976f24a144SLois Curfman McInnes snes->nwork = 0; 598c9005455SLois Curfman McInnes snes->conv_hist_size = 0; 599c9005455SLois Curfman McInnes snes->conv_act_size = 0; 600c9005455SLois Curfman McInnes snes->conv_hist = 0; 6019b94acceSBarry Smith 6029b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 6030452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 604eed86810SBarry Smith PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 6059b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6069b94acceSBarry Smith kctx->version = 2; 6079b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6089b94acceSBarry Smith this was too large for some test cases */ 6099b94acceSBarry Smith kctx->rtol_last = 0; 6109b94acceSBarry Smith kctx->rtol_max = .9; 6119b94acceSBarry Smith kctx->gamma = 1.0; 6129b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6139b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6149b94acceSBarry Smith kctx->threshold = .1; 6159b94acceSBarry Smith kctx->lresid_last = 0; 6169b94acceSBarry Smith kctx->norm_last = 0; 6179b94acceSBarry Smith 6189b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 6199b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 6205334005bSBarry Smith 6219b94acceSBarry Smith *outsnes = snes; 6229b94acceSBarry Smith return 0; 6239b94acceSBarry Smith } 6249b94acceSBarry Smith 6259b94acceSBarry Smith /* --------------------------------------------------------------- */ 6265615d1e5SSatish Balay #undef __FUNC__ 6275eea60f9SBarry Smith #define __FUNC__ "SNESSetFunction" /* ADIC Ignore */ 6289b94acceSBarry Smith /*@C 6299b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6309b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6319b94acceSBarry Smith equations. 6329b94acceSBarry Smith 6339b94acceSBarry Smith Input Parameters: 6349b94acceSBarry Smith . snes - the SNES context 6359b94acceSBarry Smith . func - function evaluation routine 6369b94acceSBarry Smith . r - vector to store function value 6372cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 6382cd2dfdcSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6399b94acceSBarry Smith 6409b94acceSBarry Smith Calling sequence of func: 641313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Vec f,void *ctx); 6429b94acceSBarry Smith 6439b94acceSBarry Smith . x - input vector 644313e4042SLois Curfman McInnes . f - function vector 6452cd2dfdcSLois Curfman McInnes . ctx - optional user-defined function context 6469b94acceSBarry Smith 6479b94acceSBarry Smith Notes: 6489b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6499b94acceSBarry Smith $ f'(x) x = -f(x), 6509b94acceSBarry Smith $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 6519b94acceSBarry Smith 6529b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6539b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6549b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 6559b94acceSBarry Smith 6569b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6579b94acceSBarry Smith 658a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6599b94acceSBarry Smith @*/ 6605334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 6619b94acceSBarry Smith { 66277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 663e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6649b94acceSBarry Smith snes->computefunction = func; 6659b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6669b94acceSBarry Smith snes->funP = ctx; 6679b94acceSBarry Smith return 0; 6689b94acceSBarry Smith } 6699b94acceSBarry Smith 6705615d1e5SSatish Balay #undef __FUNC__ 6715615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 6729b94acceSBarry Smith /*@ 6739b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6749b94acceSBarry Smith SNESSetFunction(). 6759b94acceSBarry Smith 6769b94acceSBarry Smith Input Parameters: 6779b94acceSBarry Smith . snes - the SNES context 6789b94acceSBarry Smith . x - input vector 6799b94acceSBarry Smith 6809b94acceSBarry Smith Output Parameter: 6813638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6829b94acceSBarry Smith 6831bffabb2SLois Curfman McInnes Notes: 6849b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6859b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6869b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6879b94acceSBarry Smith 6889b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6899b94acceSBarry Smith 690a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6919b94acceSBarry Smith @*/ 6929b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 6939b94acceSBarry Smith { 6949b94acceSBarry Smith int ierr; 6959b94acceSBarry Smith 69674679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 697e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6989b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 6999b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 700ae3c334cSLois Curfman McInnes snes->nfuncs++; 7019b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 7029b94acceSBarry Smith return 0; 7039b94acceSBarry Smith } 7049b94acceSBarry Smith 7055615d1e5SSatish Balay #undef __FUNC__ 7065eea60f9SBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" /* ADIC Ignore */ 7079b94acceSBarry Smith /*@C 7089b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 7099b94acceSBarry Smith unconstrained minimization. 7109b94acceSBarry Smith 7119b94acceSBarry Smith Input Parameters: 7129b94acceSBarry Smith . snes - the SNES context 7139b94acceSBarry Smith . func - function evaluation routine 714044dda88SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 715044dda88SLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 7169b94acceSBarry Smith 7179b94acceSBarry Smith Calling sequence of func: 7189b94acceSBarry Smith . func (SNES snes,Vec x,double *f,void *ctx); 7199b94acceSBarry Smith 7209b94acceSBarry Smith . x - input vector 7219b94acceSBarry Smith . f - function 722044dda88SLois Curfman McInnes . ctx - [optional] user-defined function context 7239b94acceSBarry Smith 7249b94acceSBarry Smith Notes: 7259b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7269b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7279b94acceSBarry Smith SNESSetFunction(). 7289b94acceSBarry Smith 7299b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 7309b94acceSBarry Smith 731a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 732a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 7339b94acceSBarry Smith @*/ 7349b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 7359b94acceSBarry Smith void *ctx) 7369b94acceSBarry Smith { 73777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 738e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7399b94acceSBarry Smith snes->computeumfunction = func; 7409b94acceSBarry Smith snes->umfunP = ctx; 7419b94acceSBarry Smith return 0; 7429b94acceSBarry Smith } 7439b94acceSBarry Smith 7445615d1e5SSatish Balay #undef __FUNC__ 7455615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 7469b94acceSBarry Smith /*@ 7479b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 7489b94acceSBarry Smith set with SNESSetMinimizationFunction(). 7499b94acceSBarry Smith 7509b94acceSBarry Smith Input Parameters: 7519b94acceSBarry Smith . snes - the SNES context 7529b94acceSBarry Smith . x - input vector 7539b94acceSBarry Smith 7549b94acceSBarry Smith Output Parameter: 7559b94acceSBarry Smith . y - function value 7569b94acceSBarry Smith 7579b94acceSBarry Smith Notes: 7589b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 7599b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7609b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 761a86d99e1SLois Curfman McInnes 762a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 763a86d99e1SLois Curfman McInnes 764a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 765a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 7669b94acceSBarry Smith @*/ 7679b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 7689b94acceSBarry Smith { 7699b94acceSBarry Smith int ierr; 770e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7719b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 7729b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 773ae3c334cSLois Curfman McInnes snes->nfuncs++; 7749b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7759b94acceSBarry Smith return 0; 7769b94acceSBarry Smith } 7779b94acceSBarry Smith 7785615d1e5SSatish Balay #undef __FUNC__ 7795eea60f9SBarry Smith #define __FUNC__ "SNESSetGradient" /* ADIC Ignore */ 7809b94acceSBarry Smith /*@C 7819b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 7829b94acceSBarry Smith vector for use by the SNES routines. 7839b94acceSBarry Smith 7849b94acceSBarry Smith Input Parameters: 7859b94acceSBarry Smith . snes - the SNES context 7869b94acceSBarry Smith . func - function evaluation routine 787044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 788044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 7899b94acceSBarry Smith . r - vector to store gradient value 7909b94acceSBarry Smith 7919b94acceSBarry Smith Calling sequence of func: 7929b94acceSBarry Smith . func (SNES, Vec x, Vec g, void *ctx); 7939b94acceSBarry Smith 7949b94acceSBarry Smith . x - input vector 7959b94acceSBarry Smith . g - gradient vector 796044dda88SLois Curfman McInnes . ctx - optional user-defined gradient context 7979b94acceSBarry Smith 7989b94acceSBarry Smith Notes: 7999b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 8009b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 8019b94acceSBarry Smith SNESSetFunction(). 8029b94acceSBarry Smith 8039b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 8049b94acceSBarry Smith 805a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 806a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 8079b94acceSBarry Smith @*/ 80874679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 8099b94acceSBarry Smith { 81077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 811e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8129b94acceSBarry Smith snes->computefunction = func; 8139b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 8149b94acceSBarry Smith snes->funP = ctx; 8159b94acceSBarry Smith return 0; 8169b94acceSBarry Smith } 8179b94acceSBarry Smith 8185615d1e5SSatish Balay #undef __FUNC__ 8195615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 8209b94acceSBarry Smith /*@ 821a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 822a86d99e1SLois Curfman McInnes SNESSetGradient(). 8239b94acceSBarry Smith 8249b94acceSBarry Smith Input Parameters: 8259b94acceSBarry Smith . snes - the SNES context 8269b94acceSBarry Smith . x - input vector 8279b94acceSBarry Smith 8289b94acceSBarry Smith Output Parameter: 8299b94acceSBarry Smith . y - gradient vector 8309b94acceSBarry Smith 8319b94acceSBarry Smith Notes: 8329b94acceSBarry Smith SNESComputeGradient() is valid only for 8339b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 8349b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 835a86d99e1SLois Curfman McInnes 836a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 837a86d99e1SLois Curfman McInnes 838a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 839a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 8409b94acceSBarry Smith @*/ 8419b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 8429b94acceSBarry Smith { 8439b94acceSBarry Smith int ierr; 84474679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 845e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8469b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 8479b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 8489b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 8499b94acceSBarry Smith return 0; 8509b94acceSBarry Smith } 8519b94acceSBarry Smith 8525615d1e5SSatish Balay #undef __FUNC__ 8535615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 85462fef451SLois Curfman McInnes /*@ 85562fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 85662fef451SLois Curfman McInnes set with SNESSetJacobian(). 85762fef451SLois Curfman McInnes 85862fef451SLois Curfman McInnes Input Parameters: 85962fef451SLois Curfman McInnes . snes - the SNES context 86062fef451SLois Curfman McInnes . x - input vector 86162fef451SLois Curfman McInnes 86262fef451SLois Curfman McInnes Output Parameters: 86362fef451SLois Curfman McInnes . A - Jacobian matrix 86462fef451SLois Curfman McInnes . B - optional preconditioning matrix 86562fef451SLois Curfman McInnes . flag - flag indicating matrix structure 86662fef451SLois Curfman McInnes 86762fef451SLois Curfman McInnes Notes: 86862fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 86962fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 87062fef451SLois Curfman McInnes 871dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 872dc5a77f8SLois Curfman McInnes flag parameter. 87362fef451SLois Curfman McInnes 87462fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 87562fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 87662fef451SLois Curfman McInnes methods is SNESComputeJacobian(). 87762fef451SLois Curfman McInnes 87862fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 87962fef451SLois Curfman McInnes 88062fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 88162fef451SLois Curfman McInnes @*/ 8829b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8839b94acceSBarry Smith { 8849b94acceSBarry Smith int ierr; 88574679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 886e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 8879b94acceSBarry Smith if (!snes->computejacobian) return 0; 8889b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 889c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 8909b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 8919b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 8926d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 89377c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 89477c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 8959b94acceSBarry Smith return 0; 8969b94acceSBarry Smith } 8979b94acceSBarry Smith 8985615d1e5SSatish Balay #undef __FUNC__ 8995615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 90062fef451SLois Curfman McInnes /*@ 90162fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 90262fef451SLois Curfman McInnes set with SNESSetHessian(). 90362fef451SLois Curfman McInnes 90462fef451SLois Curfman McInnes Input Parameters: 90562fef451SLois Curfman McInnes . snes - the SNES context 90662fef451SLois Curfman McInnes . x - input vector 90762fef451SLois Curfman McInnes 90862fef451SLois Curfman McInnes Output Parameters: 90962fef451SLois Curfman McInnes . A - Hessian matrix 91062fef451SLois Curfman McInnes . B - optional preconditioning matrix 91162fef451SLois Curfman McInnes . flag - flag indicating matrix structure 91262fef451SLois Curfman McInnes 91362fef451SLois Curfman McInnes Notes: 91462fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 91562fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 91662fef451SLois Curfman McInnes 917dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 918dc5a77f8SLois Curfman McInnes flag parameter. 91962fef451SLois Curfman McInnes 92062fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 92162fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 92262fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 92362fef451SLois Curfman McInnes 92462fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 92562fef451SLois Curfman McInnes 926a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 927a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 92862fef451SLois Curfman McInnes @*/ 92962fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 9309b94acceSBarry Smith { 9319b94acceSBarry Smith int ierr; 93274679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 933e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 9349b94acceSBarry Smith if (!snes->computejacobian) return 0; 93562fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 9360f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 93762fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 93862fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 9390f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 94077c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 94177c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9429b94acceSBarry Smith return 0; 9439b94acceSBarry Smith } 9449b94acceSBarry Smith 9455615d1e5SSatish Balay #undef __FUNC__ 9465eea60f9SBarry Smith #define __FUNC__ "SNESSetJacobian" /* ADIC Ignore */ 9479b94acceSBarry Smith /*@C 9489b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 949044dda88SLois Curfman McInnes location to store the matrix. 9509b94acceSBarry Smith 9519b94acceSBarry Smith Input Parameters: 9529b94acceSBarry Smith . snes - the SNES context 9539b94acceSBarry Smith . A - Jacobian matrix 9549b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 9559b94acceSBarry Smith . func - Jacobian evaluation routine 9562cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 9572cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 9589b94acceSBarry Smith 9599b94acceSBarry Smith Calling sequence of func: 960313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9619b94acceSBarry Smith 9629b94acceSBarry Smith . x - input vector 9639b94acceSBarry Smith . A - Jacobian matrix 9649b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 965ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 966ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 9672cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Jacobian context 9689b94acceSBarry Smith 9699b94acceSBarry Smith Notes: 970dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 9712cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 972ac21db08SLois Curfman McInnes 973ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9749b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9759b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9769b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9779b94acceSBarry Smith throughout the global iterations. 9789b94acceSBarry Smith 9799b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9809b94acceSBarry Smith 981ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 9829b94acceSBarry Smith @*/ 9839b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 9849b94acceSBarry Smith MatStructure*,void*),void *ctx) 9859b94acceSBarry Smith { 98677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 98774679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 988e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 9899b94acceSBarry Smith snes->computejacobian = func; 9909b94acceSBarry Smith snes->jacP = ctx; 9919b94acceSBarry Smith snes->jacobian = A; 9929b94acceSBarry Smith snes->jacobian_pre = B; 9939b94acceSBarry Smith return 0; 9949b94acceSBarry Smith } 99562fef451SLois Curfman McInnes 9965615d1e5SSatish Balay #undef __FUNC__ 9975eea60f9SBarry Smith #define __FUNC__ "SNESGetJacobian" /* ADIC Ignore */ 998b4fd4287SBarry Smith /*@ 999b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1000b4fd4287SBarry Smith provided context for evaluating the Jacobian. 1001b4fd4287SBarry Smith 1002b4fd4287SBarry Smith Input Parameter: 1003b4fd4287SBarry Smith . snes - the nonlinear solver context 1004b4fd4287SBarry Smith 1005b4fd4287SBarry Smith Output Parameters: 1006b4fd4287SBarry Smith . A - location to stash Jacobian matrix (or PETSC_NULL) 1007b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 1008b4fd4287SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 1009b4fd4287SBarry Smith 1010b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 1011b4fd4287SBarry Smith @*/ 1012b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1013b4fd4287SBarry Smith { 101474679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 1015e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 1016b4fd4287SBarry Smith if (A) *A = snes->jacobian; 1017b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 1018b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 1019b4fd4287SBarry Smith return 0; 1020b4fd4287SBarry Smith } 1021b4fd4287SBarry Smith 10225615d1e5SSatish Balay #undef __FUNC__ 10235eea60f9SBarry Smith #define __FUNC__ "SNESSetHessian" /* ADIC Ignore */ 10249b94acceSBarry Smith /*@C 10259b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 1026044dda88SLois Curfman McInnes location to store the matrix. 10279b94acceSBarry Smith 10289b94acceSBarry Smith Input Parameters: 10299b94acceSBarry Smith . snes - the SNES context 10309b94acceSBarry Smith . A - Hessian matrix 10319b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 10329b94acceSBarry Smith . func - Jacobian evaluation routine 1033313e4042SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 1034313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 10359b94acceSBarry Smith 10369b94acceSBarry Smith Calling sequence of func: 1037313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 10389b94acceSBarry Smith 10399b94acceSBarry Smith . x - input vector 10409b94acceSBarry Smith . A - Hessian matrix 10419b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1042ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1043ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 10442cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Hessian context 10459b94acceSBarry Smith 10469b94acceSBarry Smith Notes: 1047dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10482cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1049ac21db08SLois Curfman McInnes 10509b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 10519b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 10529b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10539b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10549b94acceSBarry Smith throughout the global iterations. 10559b94acceSBarry Smith 10569b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 10579b94acceSBarry Smith 1058ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 10599b94acceSBarry Smith @*/ 10609b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 10619b94acceSBarry Smith MatStructure*,void*),void *ctx) 10629b94acceSBarry Smith { 106377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 106474679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1065e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 10669b94acceSBarry Smith snes->computejacobian = func; 10679b94acceSBarry Smith snes->jacP = ctx; 10689b94acceSBarry Smith snes->jacobian = A; 10699b94acceSBarry Smith snes->jacobian_pre = B; 10709b94acceSBarry Smith return 0; 10719b94acceSBarry Smith } 10729b94acceSBarry Smith 10735615d1e5SSatish Balay #undef __FUNC__ 10745eea60f9SBarry Smith #define __FUNC__ "SNESGetHessian" /* ADIC Ignore */ 107562fef451SLois Curfman McInnes /*@ 107662fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 107762fef451SLois Curfman McInnes provided context for evaluating the Hessian. 107862fef451SLois Curfman McInnes 107962fef451SLois Curfman McInnes Input Parameter: 108062fef451SLois Curfman McInnes . snes - the nonlinear solver context 108162fef451SLois Curfman McInnes 108262fef451SLois Curfman McInnes Output Parameters: 108362fef451SLois Curfman McInnes . A - location to stash Hessian matrix (or PETSC_NULL) 108462fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 108562fef451SLois Curfman McInnes . ctx - location to stash Hessian ctx (or PETSC_NULL) 108662fef451SLois Curfman McInnes 108762fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 108862fef451SLois Curfman McInnes @*/ 108962fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 109062fef451SLois Curfman McInnes { 109174679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1092e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 109362fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 109462fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 109562fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 109662fef451SLois Curfman McInnes return 0; 109762fef451SLois Curfman McInnes } 109862fef451SLois Curfman McInnes 10999b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 11009b94acceSBarry Smith 11015615d1e5SSatish Balay #undef __FUNC__ 11025615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 11039b94acceSBarry Smith /*@ 11049b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1105272ac6f2SLois Curfman McInnes of a nonlinear solver. 11069b94acceSBarry Smith 11079b94acceSBarry Smith Input Parameter: 11089b94acceSBarry Smith . snes - the SNES context 11098ddd3da0SLois Curfman McInnes . x - the solution vector 11109b94acceSBarry Smith 1111272ac6f2SLois Curfman McInnes Notes: 1112272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1113272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1114272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1115272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1116272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1117272ac6f2SLois Curfman McInnes 11189b94acceSBarry Smith .keywords: SNES, nonlinear, setup 11199b94acceSBarry Smith 11209b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 11219b94acceSBarry Smith @*/ 11228ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 11239b94acceSBarry Smith { 1124272ac6f2SLois Curfman McInnes int ierr, flg; 112577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 112677c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 11278ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 11289b94acceSBarry Smith 1129c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1130c1f60f51SBarry Smith /* 1131c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1132dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1133c1f60f51SBarry Smith */ 1134c1f60f51SBarry Smith if (flg) { 1135c1f60f51SBarry Smith Mat J; 1136c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1137c1f60f51SBarry Smith PLogObjectParent(snes,J); 1138c1f60f51SBarry Smith snes->mfshell = J; 1139c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1140c1f60f51SBarry Smith snes->jacobian = J; 1141c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1142c1f60f51SBarry Smith } 1143c1f60f51SBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1144c1f60f51SBarry Smith snes->jacobian = J; 1145c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1146e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free operator option"); 1147c1f60f51SBarry Smith } 1148272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1149c1f60f51SBarry Smith /* 1150dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1151c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1152c1f60f51SBarry Smith */ 115331615d04SLois Curfman McInnes if (flg) { 1154272ac6f2SLois Curfman McInnes Mat J; 1155272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1156272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1157272ac6f2SLois Curfman McInnes snes->mfshell = J; 115883e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 115983e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 116094a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 116183e56afdSLois Curfman McInnes } 116283e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 116383e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 116494a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1165e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free option"); 1166272ac6f2SLois Curfman McInnes } 11679b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1168e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1169e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1170e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first"); 1171e3372554SBarry Smith if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector"); 1172a703fe33SLois Curfman McInnes 1173a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 117440191667SLois Curfman McInnes if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1175a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1176a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1177a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1178a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1179a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1180a703fe33SLois Curfman McInnes } 11819b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1182e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first"); 1183e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first"); 118437fcc0dbSBarry Smith if (!snes->computeumfunction) 1185e3372554SBarry Smith SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first"); 1186e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first"); 1187e3372554SBarry Smith } else SETERRQ(1,0,"Unknown method class"); 1188a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1189a703fe33SLois Curfman McInnes snes->setup_called = 1; 1190a703fe33SLois Curfman McInnes return 0; 11919b94acceSBarry Smith } 11929b94acceSBarry Smith 11935615d1e5SSatish Balay #undef __FUNC__ 11945eea60f9SBarry Smith #define __FUNC__ "SNESDestroy" /* ADIC Ignore */ 11959b94acceSBarry Smith /*@C 11969b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 11979b94acceSBarry Smith with SNESCreate(). 11989b94acceSBarry Smith 11999b94acceSBarry Smith Input Parameter: 12009b94acceSBarry Smith . snes - the SNES context 12019b94acceSBarry Smith 12029b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 12039b94acceSBarry Smith 120463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 12059b94acceSBarry Smith @*/ 12069b94acceSBarry Smith int SNESDestroy(SNES snes) 12079b94acceSBarry Smith { 12089b94acceSBarry Smith int ierr; 120977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12109750a799SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);} 12110452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 12129b94acceSBarry Smith if (snes->mfshell) MatDestroy(snes->mfshell); 12139b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 12143e2e452bSBarry Smith if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 12156f24a144SLois Curfman McInnes if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 12169b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 12170452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 12189b94acceSBarry Smith return 0; 12199b94acceSBarry Smith } 12209b94acceSBarry Smith 12219b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 12229b94acceSBarry Smith 12235615d1e5SSatish Balay #undef __FUNC__ 12245615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 12259b94acceSBarry Smith /*@ 1226d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 12279b94acceSBarry Smith 12289b94acceSBarry Smith Input Parameters: 12299b94acceSBarry Smith . snes - the SNES context 123033174efeSLois Curfman McInnes . atol - absolute convergence tolerance 123133174efeSLois Curfman McInnes . rtol - relative convergence tolerance 123233174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 123333174efeSLois Curfman McInnes of the change in the solution between steps 123433174efeSLois Curfman McInnes . maxit - maximum number of iterations 123533174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 12369b94acceSBarry Smith 123733174efeSLois Curfman McInnes Options Database Keys: 123833174efeSLois Curfman McInnes $ -snes_atol <atol> 123933174efeSLois Curfman McInnes $ -snes_rtol <rtol> 124033174efeSLois Curfman McInnes $ -snes_stol <stol> 124133174efeSLois Curfman McInnes $ -snes_max_it <maxit> 124233174efeSLois Curfman McInnes $ -snes_max_funcs <maxf> 12439b94acceSBarry Smith 1244d7a720efSLois Curfman McInnes Notes: 12459b94acceSBarry Smith The default maximum number of iterations is 50. 12469b94acceSBarry Smith The default maximum number of function evaluations is 1000. 12479b94acceSBarry Smith 124833174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 12499b94acceSBarry Smith 1250d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 12519b94acceSBarry Smith @*/ 1252d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 12539b94acceSBarry Smith { 125477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1255d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1256d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1257d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1258d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1259d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 12609b94acceSBarry Smith return 0; 12619b94acceSBarry Smith } 12629b94acceSBarry Smith 12635615d1e5SSatish Balay #undef __FUNC__ 12645615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 12659b94acceSBarry Smith /*@ 126633174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 126733174efeSLois Curfman McInnes 126833174efeSLois Curfman McInnes Input Parameters: 126933174efeSLois Curfman McInnes . snes - the SNES context 127033174efeSLois Curfman McInnes . atol - absolute convergence tolerance 127133174efeSLois Curfman McInnes . rtol - relative convergence tolerance 127233174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 127333174efeSLois Curfman McInnes of the change in the solution between steps 127433174efeSLois Curfman McInnes . maxit - maximum number of iterations 127533174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 127633174efeSLois Curfman McInnes 127733174efeSLois Curfman McInnes Notes: 127833174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 127933174efeSLois Curfman McInnes 128033174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 128133174efeSLois Curfman McInnes 128233174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 128333174efeSLois Curfman McInnes @*/ 128433174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 128533174efeSLois Curfman McInnes { 128633174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 128733174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 128833174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 128933174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 129033174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 129133174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 129233174efeSLois Curfman McInnes return 0; 129333174efeSLois Curfman McInnes } 129433174efeSLois Curfman McInnes 12955615d1e5SSatish Balay #undef __FUNC__ 12965615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 129733174efeSLois Curfman McInnes /*@ 12989b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 12999b94acceSBarry Smith 13009b94acceSBarry Smith Input Parameters: 13019b94acceSBarry Smith . snes - the SNES context 13029b94acceSBarry Smith . tol - tolerance 13039b94acceSBarry Smith 13049b94acceSBarry Smith Options Database Key: 1305d7a720efSLois Curfman McInnes $ -snes_trtol <tol> 13069b94acceSBarry Smith 13079b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 13089b94acceSBarry Smith 1309d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 13109b94acceSBarry Smith @*/ 13119b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 13129b94acceSBarry Smith { 131377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13149b94acceSBarry Smith snes->deltatol = tol; 13159b94acceSBarry Smith return 0; 13169b94acceSBarry Smith } 13179b94acceSBarry Smith 13185615d1e5SSatish Balay #undef __FUNC__ 13195615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 13209b94acceSBarry Smith /*@ 132177c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 13229b94acceSBarry Smith for unconstrained minimization solvers. 13239b94acceSBarry Smith 13249b94acceSBarry Smith Input Parameters: 13259b94acceSBarry Smith . snes - the SNES context 13269b94acceSBarry Smith . ftol - minimum function tolerance 13279b94acceSBarry Smith 13289b94acceSBarry Smith Options Database Key: 1329d7a720efSLois Curfman McInnes $ -snes_fmin <ftol> 13309b94acceSBarry Smith 13319b94acceSBarry Smith Note: 133277c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 13339b94acceSBarry Smith methods only. 13349b94acceSBarry Smith 13359b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 13369b94acceSBarry Smith 1337d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 13389b94acceSBarry Smith @*/ 133977c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 13409b94acceSBarry Smith { 134177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13429b94acceSBarry Smith snes->fmin = ftol; 13439b94acceSBarry Smith return 0; 13449b94acceSBarry Smith } 13459b94acceSBarry Smith 13469b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 13479b94acceSBarry Smith 13485615d1e5SSatish Balay #undef __FUNC__ 13495eea60f9SBarry Smith #define __FUNC__ "SNESSetMonitor" /* ADIC Ignore */ 13509b94acceSBarry Smith /*@C 13519b94acceSBarry Smith SNESSetMonitor - Sets the function that is to be used at every 13529b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 13539b94acceSBarry Smith progress. 13549b94acceSBarry Smith 13559b94acceSBarry Smith Input Parameters: 13569b94acceSBarry Smith . snes - the SNES context 13579b94acceSBarry Smith . func - monitoring routine 1358044dda88SLois Curfman McInnes . mctx - [optional] user-defined context for private data for the 1359044dda88SLois Curfman McInnes monitor routine (may be PETSC_NULL) 13609b94acceSBarry Smith 13619b94acceSBarry Smith Calling sequence of func: 1362682d7d0cSBarry Smith int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 13639b94acceSBarry Smith 13649b94acceSBarry Smith $ snes - the SNES context 13659b94acceSBarry Smith $ its - iteration number 1366044dda88SLois Curfman McInnes $ mctx - [optional] monitoring context 13679b94acceSBarry Smith $ 13689b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13699b94acceSBarry Smith $ norm - 2-norm function value (may be estimated) 13709b94acceSBarry Smith $ 13719b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13729b94acceSBarry Smith $ norm - 2-norm gradient value (may be estimated) 13739b94acceSBarry Smith 13749665c990SLois Curfman McInnes Options Database Keys: 13759665c990SLois Curfman McInnes $ -snes_monitor : sets SNESDefaultMonitor() 13769665c990SLois Curfman McInnes $ -snes_xmonitor : sets line graph monitor, 13779665c990SLois Curfman McInnes $ uses SNESLGMonitorCreate() 13789665c990SLois Curfman McInnes $ -snes_cancelmonitors : cancels all monitors that have 13799665c990SLois Curfman McInnes $ been hardwired into a code by 13809665c990SLois Curfman McInnes $ calls to SNESSetMonitor(), but 13819665c990SLois Curfman McInnes $ does not cancel those set via 13829665c990SLois Curfman McInnes $ the options database. 13839665c990SLois Curfman McInnes 13849665c990SLois Curfman McInnes 1385639f9d9dSBarry Smith Notes: 13866bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 13876bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 13886bc08f3fSLois Curfman McInnes order in which they were set. 1389639f9d9dSBarry Smith 13909b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 13919b94acceSBarry Smith 13929b94acceSBarry Smith .seealso: SNESDefaultMonitor() 13939b94acceSBarry Smith @*/ 139474679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 13959b94acceSBarry Smith { 1396639f9d9dSBarry Smith if (!func) { 1397639f9d9dSBarry Smith snes->numbermonitors = 0; 1398639f9d9dSBarry Smith return 0; 1399639f9d9dSBarry Smith } 1400639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1401e3372554SBarry Smith SETERRQ(1,0,"Too many monitors set"); 1402639f9d9dSBarry Smith } 1403639f9d9dSBarry Smith 1404639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1405639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 14069b94acceSBarry Smith return 0; 14079b94acceSBarry Smith } 14089b94acceSBarry Smith 14095615d1e5SSatish Balay #undef __FUNC__ 14105eea60f9SBarry Smith #define __FUNC__ "SNESSetConvergenceTest" /* ADIC Ignore */ 14119b94acceSBarry Smith /*@C 14129b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 14139b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 14149b94acceSBarry Smith 14159b94acceSBarry Smith Input Parameters: 14169b94acceSBarry Smith . snes - the SNES context 14179b94acceSBarry Smith . func - routine to test for convergence 1418044dda88SLois Curfman McInnes . cctx - [optional] context for private data for the convergence routine 1419044dda88SLois Curfman McInnes (may be PETSC_NULL) 14209b94acceSBarry Smith 14219b94acceSBarry Smith Calling sequence of func: 14229b94acceSBarry Smith int func (SNES snes,double xnorm,double gnorm, 14239b94acceSBarry Smith double f,void *cctx) 14249b94acceSBarry Smith 14259b94acceSBarry Smith $ snes - the SNES context 1426044dda88SLois Curfman McInnes $ cctx - [optional] convergence context 14279b94acceSBarry Smith $ xnorm - 2-norm of current iterate 14289b94acceSBarry Smith $ 14299b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 14309b94acceSBarry Smith $ gnorm - 2-norm of current step 14319b94acceSBarry Smith $ f - 2-norm of function 14329b94acceSBarry Smith $ 14339b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 14349b94acceSBarry Smith $ gnorm - 2-norm of current gradient 14359b94acceSBarry Smith $ f - function value 14369b94acceSBarry Smith 14379b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 14389b94acceSBarry Smith 143940191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 144040191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 14419b94acceSBarry Smith @*/ 144274679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 14439b94acceSBarry Smith { 14449b94acceSBarry Smith (snes)->converged = func; 14459b94acceSBarry Smith (snes)->cnvP = cctx; 14469b94acceSBarry Smith return 0; 14479b94acceSBarry Smith } 14489b94acceSBarry Smith 14495615d1e5SSatish Balay #undef __FUNC__ 14505eea60f9SBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" /* ADIC Ignore */ 1451c9005455SLois Curfman McInnes /*@ 1452c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1453c9005455SLois Curfman McInnes 1454c9005455SLois Curfman McInnes Input Parameters: 1455c9005455SLois Curfman McInnes . snes - iterative context obtained from SNESCreate() 1456c9005455SLois Curfman McInnes . a - array to hold history 1457c9005455SLois Curfman McInnes . na - size of a 1458c9005455SLois Curfman McInnes 1459c9005455SLois Curfman McInnes Notes: 1460c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1461c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1462c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1463c9005455SLois Curfman McInnes at each step. 1464c9005455SLois Curfman McInnes 1465c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1466c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1467c9005455SLois Curfman McInnes during the section of code that is being timed. 1468c9005455SLois Curfman McInnes 1469c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1470c9005455SLois Curfman McInnes @*/ 1471c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na) 1472c9005455SLois Curfman McInnes { 1473c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1474c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1475c9005455SLois Curfman McInnes snes->conv_hist = a; 1476c9005455SLois Curfman McInnes snes->conv_hist_size = na; 1477c9005455SLois Curfman McInnes return 0; 1478c9005455SLois Curfman McInnes } 1479c9005455SLois Curfman McInnes 1480c9005455SLois Curfman McInnes #undef __FUNC__ 14815615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 14829b94acceSBarry Smith /* 14839b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 14849b94acceSBarry Smith positive parameter delta. 14859b94acceSBarry Smith 14869b94acceSBarry Smith Input Parameters: 14879b94acceSBarry Smith . snes - the SNES context 14889b94acceSBarry Smith . y - approximate solution of linear system 14899b94acceSBarry Smith . fnorm - 2-norm of current function 14909b94acceSBarry Smith . delta - trust region size 14919b94acceSBarry Smith 14929b94acceSBarry Smith Output Parameters: 14939b94acceSBarry Smith . gpnorm - predicted function norm at the new point, assuming local 14949b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 14959b94acceSBarry Smith region, and exceeds zero otherwise. 14969b94acceSBarry Smith . ynorm - 2-norm of the step 14979b94acceSBarry Smith 14989b94acceSBarry Smith Note: 149940191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 15009b94acceSBarry Smith is set to be the maximum allowable step size. 15019b94acceSBarry Smith 15029b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 15039b94acceSBarry Smith */ 15049b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 15059b94acceSBarry Smith double *gpnorm,double *ynorm) 15069b94acceSBarry Smith { 15079b94acceSBarry Smith double norm; 15089b94acceSBarry Smith Scalar cnorm; 1509cddf8d76SBarry Smith VecNorm(y,NORM_2, &norm ); 15109b94acceSBarry Smith if (norm > *delta) { 15119b94acceSBarry Smith norm = *delta/norm; 15129b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 15139b94acceSBarry Smith cnorm = norm; 15149b94acceSBarry Smith VecScale( &cnorm, y ); 15159b94acceSBarry Smith *ynorm = *delta; 15169b94acceSBarry Smith } else { 15179b94acceSBarry Smith *gpnorm = 0.0; 15189b94acceSBarry Smith *ynorm = norm; 15199b94acceSBarry Smith } 15209b94acceSBarry Smith return 0; 15219b94acceSBarry Smith } 15229b94acceSBarry Smith 15235615d1e5SSatish Balay #undef __FUNC__ 15245615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 15259b94acceSBarry Smith /*@ 15269b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 152763a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 15289b94acceSBarry Smith 15299b94acceSBarry Smith Input Parameter: 15309b94acceSBarry Smith . snes - the SNES context 15318ddd3da0SLois Curfman McInnes . x - the solution vector 15329b94acceSBarry Smith 15339b94acceSBarry Smith Output Parameter: 15349b94acceSBarry Smith its - number of iterations until termination 15359b94acceSBarry Smith 15368ddd3da0SLois Curfman McInnes Note: 15378ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 15388ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 15398ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 15408ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 15418ddd3da0SLois Curfman McInnes 15429b94acceSBarry Smith .keywords: SNES, nonlinear, solve 15439b94acceSBarry Smith 154463a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 15459b94acceSBarry Smith @*/ 15468ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 15479b94acceSBarry Smith { 15483c7409f5SSatish Balay int ierr, flg; 1549052efed2SBarry Smith 155077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 155174679c65SBarry Smith PetscValidIntPointer(its); 1552c4fc05e7SBarry Smith if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1553c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 15549b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 1555c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 15569b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 15579b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 15583c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 15596d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 15609b94acceSBarry Smith return 0; 15619b94acceSBarry Smith } 15629b94acceSBarry Smith 15639b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 156437fcc0dbSBarry Smith static NRList *__SNESList = 0; 15659b94acceSBarry Smith 15665615d1e5SSatish Balay #undef __FUNC__ 15675615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 15689b94acceSBarry Smith /*@ 15694b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 15709b94acceSBarry Smith 15719b94acceSBarry Smith Input Parameters: 15729b94acceSBarry Smith . snes - the SNES context 15739b94acceSBarry Smith . method - a known method 15749b94acceSBarry Smith 1575ae12b187SLois Curfman McInnes Options Database Command: 1576ae12b187SLois Curfman McInnes $ -snes_type <method> 1577ae12b187SLois Curfman McInnes $ Use -help for a list of available methods 1578ae12b187SLois Curfman McInnes $ (for instance, ls or tr) 1579ae12b187SLois Curfman McInnes 15809b94acceSBarry Smith Notes: 15819b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 15829b94acceSBarry Smith $ Systems of nonlinear equations: 158340191667SLois Curfman McInnes $ SNES_EQ_LS - Newton's method with line search 158440191667SLois Curfman McInnes $ SNES_EQ_TR - Newton's method with trust region 15859b94acceSBarry Smith $ Unconstrained minimization: 158640191667SLois Curfman McInnes $ SNES_UM_TR - Newton's method with trust region 158740191667SLois Curfman McInnes $ SNES_UM_LS - Newton's method with line search 15889b94acceSBarry Smith 1589ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1590ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1591ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1592ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1593ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1594ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1595ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1596ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1597ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1598ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1599a703fe33SLois Curfman McInnes 1600f536c699SSatish Balay .keywords: SNES, set, method 16019b94acceSBarry Smith @*/ 16024b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 16039b94acceSBarry Smith { 160484cb2905SBarry Smith int ierr; 16059b94acceSBarry Smith int (*r)(SNES); 160677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16077d1a2b2bSBarry Smith if (snes->type == (int) method) return 0; 16087d1a2b2bSBarry Smith 16099b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 161084cb2905SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1611e3372554SBarry Smith if (!__SNESList) {SETERRQ(1,0,"Could not get methods");} 161237fcc0dbSBarry Smith r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1613e3372554SBarry Smith if (!r) {SETERRQ(1,0,"Unknown method");} 16140452661fSBarry Smith if (snes->data) PetscFree(snes->data); 16159b94acceSBarry Smith snes->set_method_called = 1; 16169b94acceSBarry Smith return (*r)(snes); 16179b94acceSBarry Smith } 16189b94acceSBarry Smith 16199b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16205615d1e5SSatish Balay #undef __FUNC__ 16215eea60f9SBarry Smith #define __FUNC__ "SNESRegister" /* ADIC Ignore */ 16229b94acceSBarry Smith /*@C 16239b94acceSBarry Smith SNESRegister - Adds the method to the nonlinear solver package, given 16244b0e389bSBarry Smith a function pointer and a nonlinear solver name of the type SNESType. 16259b94acceSBarry Smith 16269b94acceSBarry Smith Input Parameters: 16272d872ea7SLois Curfman McInnes . name - either a predefined name such as SNES_EQ_LS, or SNES_NEW 16282d872ea7SLois Curfman McInnes to indicate a new user-defined solver 162940191667SLois Curfman McInnes . sname - corresponding string for name 16309b94acceSBarry Smith . create - routine to create method context 16319b94acceSBarry Smith 163284cb2905SBarry Smith Output Parameter: 163384cb2905SBarry Smith . oname - type associated with this new solver 163484cb2905SBarry Smith 16352d872ea7SLois Curfman McInnes Notes: 16362d872ea7SLois Curfman McInnes Multiple user-defined nonlinear solvers can be added by calling 16372d872ea7SLois Curfman McInnes SNESRegister() with the input parameter "name" set to be SNES_NEW; 16382d872ea7SLois Curfman McInnes each call will return a unique solver type in the output 16392d872ea7SLois Curfman McInnes parameter "oname". 16402d872ea7SLois Curfman McInnes 16419b94acceSBarry Smith .keywords: SNES, nonlinear, register 16429b94acceSBarry Smith 16439b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy() 16449b94acceSBarry Smith @*/ 164584cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES)) 16469b94acceSBarry Smith { 16479b94acceSBarry Smith int ierr; 164884cb2905SBarry Smith static int numberregistered = 0; 164984cb2905SBarry Smith 1650d252947aSBarry Smith if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++); 165184cb2905SBarry Smith 165284cb2905SBarry Smith if (oname) *oname = name; 165337fcc0dbSBarry Smith if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 165484cb2905SBarry Smith NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create ); 16559b94acceSBarry Smith return 0; 16569b94acceSBarry Smith } 1657a847f771SSatish Balay 16589b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16595615d1e5SSatish Balay #undef __FUNC__ 16605eea60f9SBarry Smith #define __FUNC__ "SNESRegisterDestroy" /* ADIC Ignore */ 16619b94acceSBarry Smith /*@C 16629b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 16639b94acceSBarry Smith registered by SNESRegister(). 16649b94acceSBarry Smith 16659b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 16669b94acceSBarry Smith 16679b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 16689b94acceSBarry Smith @*/ 16699b94acceSBarry Smith int SNESRegisterDestroy() 16709b94acceSBarry Smith { 167137fcc0dbSBarry Smith if (__SNESList) { 167237fcc0dbSBarry Smith NRDestroy( __SNESList ); 167337fcc0dbSBarry Smith __SNESList = 0; 16749b94acceSBarry Smith } 167584cb2905SBarry Smith SNESRegisterAllCalled = 0; 16769b94acceSBarry Smith return 0; 16779b94acceSBarry Smith } 16789b94acceSBarry Smith 16795615d1e5SSatish Balay #undef __FUNC__ 16805eea60f9SBarry Smith #define __FUNC__ "SNESGetTypeFromOptions_Private" /* ADIC Ignore */ 16819b94acceSBarry Smith /* 16824b0e389bSBarry Smith SNESGetTypeFromOptions_Private - Sets the selected method from the 16839b94acceSBarry Smith options database. 16849b94acceSBarry Smith 16859b94acceSBarry Smith Input Parameter: 16869b94acceSBarry Smith . ctx - the SNES context 16879b94acceSBarry Smith 16889b94acceSBarry Smith Output Parameter: 16899b94acceSBarry Smith . method - solver method 16909b94acceSBarry Smith 16919b94acceSBarry Smith Returns: 16929b94acceSBarry Smith Returns 1 if the method is found; 0 otherwise. 16939b94acceSBarry Smith 16949b94acceSBarry Smith Options Database Key: 16954b0e389bSBarry Smith $ -snes_type method 16969b94acceSBarry Smith */ 1697052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 16989b94acceSBarry Smith { 1699052efed2SBarry Smith int ierr; 17009b94acceSBarry Smith char sbuf[50]; 17015baf8537SBarry Smith 1702052efed2SBarry Smith ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1703052efed2SBarry Smith if (*flg) { 1704052efed2SBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 170537fcc0dbSBarry Smith *method = (SNESType)NRFindID( __SNESList, sbuf ); 17069b94acceSBarry Smith } 17079b94acceSBarry Smith return 0; 17089b94acceSBarry Smith } 17099b94acceSBarry Smith 17105615d1e5SSatish Balay #undef __FUNC__ 17115eea60f9SBarry Smith #define __FUNC__ "SNESGetType" /* ADIC Ignore */ 17129b94acceSBarry Smith /*@C 17139a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 17149b94acceSBarry Smith 17159b94acceSBarry Smith Input Parameter: 17164b0e389bSBarry Smith . snes - nonlinear solver context 17179b94acceSBarry Smith 17189b94acceSBarry Smith Output Parameter: 17199a28b0a6SLois Curfman McInnes . method - SNES method (or use PETSC_NULL) 17209a28b0a6SLois Curfman McInnes . name - name of SNES method (or use PETSC_NULL) 17219b94acceSBarry Smith 17229b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 17239b94acceSBarry Smith @*/ 17244b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name) 17259b94acceSBarry Smith { 172606a9b0adSLois Curfman McInnes int ierr; 172737fcc0dbSBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 17284b0e389bSBarry Smith if (method) *method = (SNESType) snes->type; 172937fcc0dbSBarry Smith if (name) *name = NRFindName( __SNESList, (int) snes->type ); 17309b94acceSBarry Smith return 0; 17319b94acceSBarry Smith } 17329b94acceSBarry Smith 17339b94acceSBarry Smith #include <stdio.h> 17345615d1e5SSatish Balay #undef __FUNC__ 17355eea60f9SBarry Smith #define __FUNC__ "SNESPrintTypes_Private" /* ADIC Ignore */ 17369b94acceSBarry Smith /* 17374b0e389bSBarry Smith SNESPrintTypes_Private - Prints the SNES methods available from the 17389b94acceSBarry Smith options database. 17399b94acceSBarry Smith 17409b94acceSBarry Smith Input Parameters: 174133455573SSatish Balay . comm - communicator (usually MPI_COMM_WORLD) 17429b94acceSBarry Smith . prefix - prefix (usually "-") 17434b0e389bSBarry Smith . name - the options database name (by default "snes_type") 17449b94acceSBarry Smith */ 174533455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 17469b94acceSBarry Smith { 17479b94acceSBarry Smith FuncList *entry; 174837fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 174937fcc0dbSBarry Smith entry = __SNESList->head; 175077c4ece6SBarry Smith PetscPrintf(comm," %s%s (one of)",prefix,name); 17519b94acceSBarry Smith while (entry) { 175277c4ece6SBarry Smith PetscPrintf(comm," %s",entry->name); 17539b94acceSBarry Smith entry = entry->next; 17549b94acceSBarry Smith } 175577c4ece6SBarry Smith PetscPrintf(comm,"\n"); 17569b94acceSBarry Smith return 0; 17579b94acceSBarry Smith } 17589b94acceSBarry Smith 17595615d1e5SSatish Balay #undef __FUNC__ 17605eea60f9SBarry Smith #define __FUNC__ "SNESGetSolution" /* ADIC Ignore */ 17619b94acceSBarry Smith /*@C 17629b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 17639b94acceSBarry Smith stored. 17649b94acceSBarry Smith 17659b94acceSBarry Smith Input Parameter: 17669b94acceSBarry Smith . snes - the SNES context 17679b94acceSBarry Smith 17689b94acceSBarry Smith Output Parameter: 17699b94acceSBarry Smith . x - the solution 17709b94acceSBarry Smith 17719b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 17729b94acceSBarry Smith 1773a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 17749b94acceSBarry Smith @*/ 17759b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 17769b94acceSBarry Smith { 177777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17789b94acceSBarry Smith *x = snes->vec_sol_always; 17799b94acceSBarry Smith return 0; 17809b94acceSBarry Smith } 17819b94acceSBarry Smith 17825615d1e5SSatish Balay #undef __FUNC__ 17835eea60f9SBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" /* ADIC Ignore */ 17849b94acceSBarry Smith /*@C 17859b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 17869b94acceSBarry Smith stored. 17879b94acceSBarry Smith 17889b94acceSBarry Smith Input Parameter: 17899b94acceSBarry Smith . snes - the SNES context 17909b94acceSBarry Smith 17919b94acceSBarry Smith Output Parameter: 17929b94acceSBarry Smith . x - the solution update 17939b94acceSBarry Smith 17949b94acceSBarry Smith Notes: 17959b94acceSBarry Smith This vector is implementation dependent. 17969b94acceSBarry Smith 17979b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 17989b94acceSBarry Smith 17999b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 18009b94acceSBarry Smith @*/ 18019b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 18029b94acceSBarry Smith { 180377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18049b94acceSBarry Smith *x = snes->vec_sol_update_always; 18059b94acceSBarry Smith return 0; 18069b94acceSBarry Smith } 18079b94acceSBarry Smith 18085615d1e5SSatish Balay #undef __FUNC__ 18095eea60f9SBarry Smith #define __FUNC__ "SNESGetFunction" /* ADIC Ignore */ 18109b94acceSBarry Smith /*@C 18113638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 18129b94acceSBarry Smith 18139b94acceSBarry Smith Input Parameter: 18149b94acceSBarry Smith . snes - the SNES context 18159b94acceSBarry Smith 18169b94acceSBarry Smith Output Parameter: 18173638b69dSLois Curfman McInnes . r - the function 18189b94acceSBarry Smith 18199b94acceSBarry Smith Notes: 18209b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 18219b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 18229b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 18239b94acceSBarry Smith 1824a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 18259b94acceSBarry Smith 182631615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 182731615d04SLois Curfman McInnes SNESGetGradient() 18289b94acceSBarry Smith @*/ 18299b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 18309b94acceSBarry Smith { 183177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1832e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0, 18337c792091SSatish Balay "For SNES_NONLINEAR_EQUATIONS only"); 18349b94acceSBarry Smith *r = snes->vec_func_always; 18359b94acceSBarry Smith return 0; 18369b94acceSBarry Smith } 18379b94acceSBarry Smith 18385615d1e5SSatish Balay #undef __FUNC__ 18395eea60f9SBarry Smith #define __FUNC__ "SNESGetGradient" /* ADIC Ignore */ 18409b94acceSBarry Smith /*@C 18413638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 18429b94acceSBarry Smith 18439b94acceSBarry Smith Input Parameter: 18449b94acceSBarry Smith . snes - the SNES context 18459b94acceSBarry Smith 18469b94acceSBarry Smith Output Parameter: 18479b94acceSBarry Smith . r - the gradient 18489b94acceSBarry Smith 18499b94acceSBarry Smith Notes: 18509b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 18519b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18529b94acceSBarry Smith SNESGetFunction(). 18539b94acceSBarry Smith 18549b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 18559b94acceSBarry Smith 185631615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 18579b94acceSBarry Smith @*/ 18589b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 18599b94acceSBarry Smith { 186077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1861e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 186263c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 18639b94acceSBarry Smith *r = snes->vec_func_always; 18649b94acceSBarry Smith return 0; 18659b94acceSBarry Smith } 18669b94acceSBarry Smith 18675615d1e5SSatish Balay #undef __FUNC__ 18685eea60f9SBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" /* ADIC Ignore */ 18699b94acceSBarry Smith /*@ 18709b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 18719b94acceSBarry Smith unconstrained minimization problems. 18729b94acceSBarry Smith 18739b94acceSBarry Smith Input Parameter: 18749b94acceSBarry Smith . snes - the SNES context 18759b94acceSBarry Smith 18769b94acceSBarry Smith Output Parameter: 18779b94acceSBarry Smith . r - the function 18789b94acceSBarry Smith 18799b94acceSBarry Smith Notes: 18809b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 18819b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18829b94acceSBarry Smith SNESGetFunction(). 18839b94acceSBarry Smith 18849b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 18859b94acceSBarry Smith 188631615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 18879b94acceSBarry Smith @*/ 18889b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 18899b94acceSBarry Smith { 189077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 189174679c65SBarry Smith PetscValidScalarPointer(r); 1892e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 189363c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 18949b94acceSBarry Smith *r = snes->fc; 18959b94acceSBarry Smith return 0; 18969b94acceSBarry Smith } 18979b94acceSBarry Smith 18985615d1e5SSatish Balay #undef __FUNC__ 18995eea60f9SBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" /* ADIC Ignore */ 19003c7409f5SSatish Balay /*@C 19013c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1902d850072dSLois Curfman McInnes SNES options in the database. 19033c7409f5SSatish Balay 19043c7409f5SSatish Balay Input Parameter: 19053c7409f5SSatish Balay . snes - the SNES context 19063c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 19073c7409f5SSatish Balay 1908d850072dSLois Curfman McInnes Notes: 1909a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1910a83b1b31SSatish Balay The first character of all runtime options is AUTOMATICALLY the 1911a83b1b31SSatish Balay hyphen. 1912d850072dSLois Curfman McInnes 19133c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1914a86d99e1SLois Curfman McInnes 1915a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19163c7409f5SSatish Balay @*/ 19173c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 19183c7409f5SSatish Balay { 19193c7409f5SSatish Balay int ierr; 19203c7409f5SSatish Balay 192177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1922639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19233c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19243c7409f5SSatish Balay return 0; 19253c7409f5SSatish Balay } 19263c7409f5SSatish Balay 19275615d1e5SSatish Balay #undef __FUNC__ 19285eea60f9SBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" /* ADIC Ignore */ 19293c7409f5SSatish Balay /*@C 1930f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1931d850072dSLois Curfman McInnes SNES options in the database. 19323c7409f5SSatish Balay 19333c7409f5SSatish Balay Input Parameter: 19343c7409f5SSatish Balay . snes - the SNES context 19353c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 19363c7409f5SSatish Balay 1937d850072dSLois Curfman McInnes Notes: 1938a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1939a83b1b31SSatish Balay The first character of all runtime options is AUTOMATICALLY the 1940a83b1b31SSatish Balay hyphen. 1941d850072dSLois Curfman McInnes 19423c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1943a86d99e1SLois Curfman McInnes 1944a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 19453c7409f5SSatish Balay @*/ 19463c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 19473c7409f5SSatish Balay { 19483c7409f5SSatish Balay int ierr; 19493c7409f5SSatish Balay 195077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1951639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19523c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19533c7409f5SSatish Balay return 0; 19543c7409f5SSatish Balay } 19553c7409f5SSatish Balay 19565615d1e5SSatish Balay #undef __FUNC__ 19575eea60f9SBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" /* ADIC Ignore */ 1958c83590e2SSatish Balay /*@ 19593c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 19603c7409f5SSatish Balay SNES options in the database. 19613c7409f5SSatish Balay 19623c7409f5SSatish Balay Input Parameter: 19633c7409f5SSatish Balay . snes - the SNES context 19643c7409f5SSatish Balay 19653c7409f5SSatish Balay Output Parameter: 19663c7409f5SSatish Balay . prefix - pointer to the prefix string used 19673c7409f5SSatish Balay 19683c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 1969a86d99e1SLois Curfman McInnes 1970a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 19713c7409f5SSatish Balay @*/ 19723c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 19733c7409f5SSatish Balay { 19743c7409f5SSatish Balay int ierr; 19753c7409f5SSatish Balay 197677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1977639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19783c7409f5SSatish Balay return 0; 19793c7409f5SSatish Balay } 19803c7409f5SSatish Balay 19813c7409f5SSatish Balay 19829b94acceSBarry Smith 19839b94acceSBarry Smith 19849b94acceSBarry Smith 1985