1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*83e2fdc7SBarry Smith static char vcid[] = "$Id: snes.c,v 1.125 1997/07/09 20:59:37 balay 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 132*83e2fdc7SBarry Smith Notes: To see all options, run your program with the -help option; 133*83e2fdc7SBarry Smith or consult the users manual. 134*83e2fdc7SBarry Smith 1359b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1369b94acceSBarry Smith 137a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1389b94acceSBarry Smith @*/ 1399b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1409b94acceSBarry Smith { 1414b0e389bSBarry Smith SNESType method; 1429b94acceSBarry Smith double tmp; 1439b94acceSBarry Smith SLES sles; 14481f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1453c7409f5SSatish Balay int version = PETSC_DEFAULT; 1469b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1479b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1489b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1499b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1509b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1519b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1529b94acceSBarry Smith 15381f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 15481f57ec7SBarry Smith 15581f57ec7SBarry Smith 15677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 157e3372554SBarry Smith if (snes->setup_called) SETERRQ(1,0,"Must call prior to SNESSetUp"); 158052efed2SBarry Smith ierr = SNESGetTypeFromOptions_Private(snes,&method,&flg); CHKERRQ(ierr); 159052efed2SBarry Smith if (flg) { 160052efed2SBarry Smith ierr = SNESSetType(snes,method); CHKERRQ(ierr); 1619b94acceSBarry Smith } 1625334005bSBarry Smith else if (!snes->set_method_called) { 1635334005bSBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 16440191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 1655334005bSBarry Smith } 1665334005bSBarry Smith else { 16740191667SLois Curfman McInnes ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 1685334005bSBarry Smith } 1695334005bSBarry Smith } 1705334005bSBarry Smith 1713c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 1723c7409f5SSatish Balay if (flg) { SNESPrintHelp(snes); } 1733c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 174d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);} 1753c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 176d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);} 1773c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 178d7a720efSLois Curfman McInnes if (flg) {SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT); } 1793c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 1803c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 181d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 182d7a720efSLois Curfman McInnes if (flg) { SNESSetTrustRegionTolerance(snes,tmp); } 183d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 184d7a720efSLois Curfman McInnes if (flg) { SNESSetMinimizationFunctionTolerance(snes,tmp); } 1853c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 1863c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 187b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 188b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 189b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 190b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 191b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 192b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 193b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 19493c39befSBarry Smith 19593c39befSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr); 19693c39befSBarry Smith if (flg) {snes->converged = 0;} 19793c39befSBarry Smith 1989b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 1999b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 2005bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 2015bbad29bSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,0,0);CHKERRQ(ierr);} 2023c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 203639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 2043c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 205d31a3109SLois Curfman McInnes if (flg) { 206639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 207d31a3109SLois Curfman McInnes } 20881f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2093c7409f5SSatish Balay if (flg) { 21017699dbbSLois Curfman McInnes int rank = 0; 211d7e8b826SBarry Smith DrawLG lg; 21217699dbbSLois Curfman McInnes MPI_Initialized(&rank); 21317699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 21417699dbbSLois Curfman McInnes if (!rank) { 21581f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2169b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 217f8c826e1SBarry Smith snes->xmonitor = lg; 2189b94acceSBarry Smith PLogObjectParent(snes,lg); 2199b94acceSBarry Smith } 2209b94acceSBarry Smith } 2213c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2223c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2239b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2249b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 22594a424c1SBarry Smith PLogInfo(snes, 22631615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 22731615d04SLois Curfman McInnes } 22831615d04SLois Curfman McInnes else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 22931615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 23031615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 23194a424c1SBarry Smith PLogInfo(snes, 23231615d04SLois Curfman McInnes "SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2339b94acceSBarry Smith } 234639f9d9dSBarry Smith 235639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 236639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 237639f9d9dSBarry Smith } 238639f9d9dSBarry Smith 2399b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2409b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 2419b94acceSBarry Smith if (!snes->setfromoptions) return 0; 2429b94acceSBarry Smith return (*snes->setfromoptions)(snes); 2439b94acceSBarry Smith } 2449b94acceSBarry Smith 2455615d1e5SSatish Balay #undef __FUNC__ 2465eea60f9SBarry Smith #define __FUNC__ "SNESPrintHelp" /* ADIC Ignore */ 2479b94acceSBarry Smith /*@ 2489b94acceSBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2499b94acceSBarry Smith 2509b94acceSBarry Smith Input Parameter: 2519b94acceSBarry Smith . snes - the SNES context 2529b94acceSBarry Smith 253a703fe33SLois Curfman McInnes Options Database Keys: 254a703fe33SLois Curfman McInnes $ -help, -h 255a703fe33SLois Curfman McInnes 2569b94acceSBarry Smith .keywords: SNES, nonlinear, help 2579b94acceSBarry Smith 258682d7d0cSBarry Smith .seealso: SNESSetFromOptions() 2599b94acceSBarry Smith @*/ 2609b94acceSBarry Smith int SNESPrintHelp(SNES snes) 2619b94acceSBarry Smith { 2626daaf66cSBarry Smith char p[64]; 2636daaf66cSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2646daaf66cSBarry Smith 26577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2666daaf66cSBarry Smith 2676daaf66cSBarry Smith PetscStrcpy(p,"-"); 2686daaf66cSBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2696daaf66cSBarry Smith 2706daaf66cSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2716daaf66cSBarry Smith 272d31a3109SLois Curfman McInnes PetscPrintf(snes->comm,"SNES options ------------------------------------------------\n"); 27333455573SSatish Balay SNESPrintTypes_Private(snes->comm,p,"snes_type"); 27477c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 275b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 276b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 277b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 278b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 279b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 280b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 2815bbad29bSBarry Smith PetscPrintf(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 2825bbad29bSBarry Smith PetscPrintf(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 283*83e2fdc7SBarry Smith PetscPrintf(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 284*83e2fdc7SBarry Smith residual norm at each iteration.\n",p); 285*83e2fdc7SBarry Smith PetscPrintf(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 286*83e2fdc7SBarry Smith residual norm for small residual norms. This is useful to conceal\n\ 287*83e2fdc7SBarry Smith meaningless digits that may be different on different machines.\n",p); 288d31a3109SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 289b18e04deSLois Curfman McInnes if (snes->type == SNES_NONLINEAR_EQUATIONS) { 29077c4ece6SBarry Smith PetscPrintf(snes->comm, 291d31a3109SLois Curfman McInnes " Options for solving systems of nonlinear equations only:\n"); 29277c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 29377c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 294ef1dfb11SLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 29593c39befSBarry Smith PetscPrintf(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p); 29677c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 29777c4ece6SBarry Smith PetscPrintf(snes->comm, 298b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 29977c4ece6SBarry Smith PetscPrintf(snes->comm, 300b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 30177c4ece6SBarry Smith PetscPrintf(snes->comm, 302b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 30377c4ece6SBarry Smith PetscPrintf(snes->comm, 304b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 30577c4ece6SBarry Smith PetscPrintf(snes->comm, 306b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 30777c4ece6SBarry Smith PetscPrintf(snes->comm, 308b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 30977c4ece6SBarry Smith PetscPrintf(snes->comm, 310b18e04deSLois Curfman McInnes " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 311b18e04deSLois Curfman McInnes } 312b18e04deSLois Curfman McInnes else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 31377c4ece6SBarry Smith PetscPrintf(snes->comm, 314d31a3109SLois Curfman McInnes " Options for solving unconstrained minimization problems only:\n"); 315b18e04deSLois Curfman McInnes PetscPrintf(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 31677c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 31777c4ece6SBarry Smith PetscPrintf(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 318b18e04deSLois Curfman McInnes } 3194537a946SLois Curfman McInnes PetscPrintf(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 32077c4ece6SBarry Smith PetscPrintf(snes->comm,"a particular method\n"); 3216daaf66cSBarry Smith if (snes->printhelp) (*snes->printhelp)(snes,p); 3229b94acceSBarry Smith return 0; 3239b94acceSBarry Smith } 324a847f771SSatish Balay 3255615d1e5SSatish Balay #undef __FUNC__ 3265eea60f9SBarry Smith #define __FUNC__ "SNESSetApplicationContext" /* ADIC Ignore */ 3279b94acceSBarry Smith /*@ 3289b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 3299b94acceSBarry Smith the nonlinear solvers. 3309b94acceSBarry Smith 3319b94acceSBarry Smith Input Parameters: 3329b94acceSBarry Smith . snes - the SNES context 3339b94acceSBarry Smith . usrP - optional user context 3349b94acceSBarry Smith 3359b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3369b94acceSBarry Smith 3379b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3389b94acceSBarry Smith @*/ 3399b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3409b94acceSBarry Smith { 34177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3429b94acceSBarry Smith snes->user = usrP; 3439b94acceSBarry Smith return 0; 3449b94acceSBarry Smith } 34574679c65SBarry Smith 3465615d1e5SSatish Balay #undef __FUNC__ 3475eea60f9SBarry Smith #define __FUNC__ "SNESGetApplicationContext" /* ADIC Ignore */ 3489b94acceSBarry Smith /*@C 3499b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3509b94acceSBarry Smith nonlinear solvers. 3519b94acceSBarry Smith 3529b94acceSBarry Smith Input Parameter: 3539b94acceSBarry Smith . snes - SNES context 3549b94acceSBarry Smith 3559b94acceSBarry Smith Output Parameter: 3569b94acceSBarry Smith . usrP - user context 3579b94acceSBarry Smith 3589b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3599b94acceSBarry Smith 3609b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3619b94acceSBarry Smith @*/ 3629b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3639b94acceSBarry Smith { 36477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3659b94acceSBarry Smith *usrP = snes->user; 3669b94acceSBarry Smith return 0; 3679b94acceSBarry Smith } 36874679c65SBarry Smith 3695615d1e5SSatish Balay #undef __FUNC__ 3705eea60f9SBarry Smith #define __FUNC__ "SNESGetIterationNumber" /* ADIC Ignore */ 3719b94acceSBarry Smith /*@ 3729b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3739b94acceSBarry Smith nonlinear solver. 3749b94acceSBarry Smith 3759b94acceSBarry Smith Input Parameter: 3769b94acceSBarry Smith . snes - SNES context 3779b94acceSBarry Smith 3789b94acceSBarry Smith Output Parameter: 3799b94acceSBarry Smith . iter - iteration number 3809b94acceSBarry Smith 3819b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3829b94acceSBarry Smith @*/ 3839b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3849b94acceSBarry Smith { 38577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 38674679c65SBarry Smith PetscValidIntPointer(iter); 3879b94acceSBarry Smith *iter = snes->iter; 3889b94acceSBarry Smith return 0; 3899b94acceSBarry Smith } 39074679c65SBarry Smith 3915615d1e5SSatish Balay #undef __FUNC__ 3925615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 3939b94acceSBarry Smith /*@ 3949b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3959b94acceSBarry Smith with SNESSSetFunction(). 3969b94acceSBarry Smith 3979b94acceSBarry Smith Input Parameter: 3989b94acceSBarry Smith . snes - SNES context 3999b94acceSBarry Smith 4009b94acceSBarry Smith Output Parameter: 4019b94acceSBarry Smith . fnorm - 2-norm of function 4029b94acceSBarry Smith 4039b94acceSBarry Smith Note: 4049b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 405a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 406a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 4079b94acceSBarry Smith 4089b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 409a86d99e1SLois Curfman McInnes 410a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 4119b94acceSBarry Smith @*/ 4129b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 4139b94acceSBarry Smith { 41477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 41574679c65SBarry Smith PetscValidScalarPointer(fnorm); 41674679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 417d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 41874679c65SBarry Smith } 4199b94acceSBarry Smith *fnorm = snes->norm; 4209b94acceSBarry Smith return 0; 4219b94acceSBarry Smith } 42274679c65SBarry Smith 4235615d1e5SSatish Balay #undef __FUNC__ 4245615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 4259b94acceSBarry Smith /*@ 4269b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 4279b94acceSBarry Smith with SNESSSetGradient(). 4289b94acceSBarry Smith 4299b94acceSBarry Smith Input Parameter: 4309b94acceSBarry Smith . snes - SNES context 4319b94acceSBarry Smith 4329b94acceSBarry Smith Output Parameter: 4339b94acceSBarry Smith . fnorm - 2-norm of gradient 4349b94acceSBarry Smith 4359b94acceSBarry Smith Note: 4369b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 437a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 438a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4399b94acceSBarry Smith 4409b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 441a86d99e1SLois Curfman McInnes 442a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4439b94acceSBarry Smith @*/ 4449b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4459b94acceSBarry Smith { 44677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 44774679c65SBarry Smith PetscValidScalarPointer(gnorm); 44874679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 449d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 45074679c65SBarry Smith } 4519b94acceSBarry Smith *gnorm = snes->norm; 4529b94acceSBarry Smith return 0; 4539b94acceSBarry Smith } 45474679c65SBarry Smith 4555615d1e5SSatish Balay #undef __FUNC__ 4565eea60f9SBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" /* ADIC Ignore */ 4579b94acceSBarry Smith /*@ 4589b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4599b94acceSBarry Smith attempted by the nonlinear solver. 4609b94acceSBarry Smith 4619b94acceSBarry Smith Input Parameter: 4629b94acceSBarry Smith . snes - SNES context 4639b94acceSBarry Smith 4649b94acceSBarry Smith Output Parameter: 4659b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4669b94acceSBarry Smith 467c96a6f78SLois Curfman McInnes Notes: 468c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 469c96a6f78SLois Curfman McInnes 4709b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4719b94acceSBarry Smith @*/ 4729b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4739b94acceSBarry Smith { 47477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 47574679c65SBarry Smith PetscValidIntPointer(nfails); 4769b94acceSBarry Smith *nfails = snes->nfailures; 4779b94acceSBarry Smith return 0; 4789b94acceSBarry Smith } 479a847f771SSatish Balay 4805615d1e5SSatish Balay #undef __FUNC__ 4815eea60f9SBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" /* ADIC Ignore */ 482c96a6f78SLois Curfman McInnes /*@ 483c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 484c96a6f78SLois Curfman McInnes used by the nonlinear solver. 485c96a6f78SLois Curfman McInnes 486c96a6f78SLois Curfman McInnes Input Parameter: 487c96a6f78SLois Curfman McInnes . snes - SNES context 488c96a6f78SLois Curfman McInnes 489c96a6f78SLois Curfman McInnes Output Parameter: 490c96a6f78SLois Curfman McInnes . lits - number of linear iterations 491c96a6f78SLois Curfman McInnes 492c96a6f78SLois Curfman McInnes Notes: 493c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 494c96a6f78SLois Curfman McInnes 495c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 496c96a6f78SLois Curfman McInnes @*/ 49786bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 498c96a6f78SLois Curfman McInnes { 499c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 500c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 501c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 502c96a6f78SLois Curfman McInnes return 0; 503c96a6f78SLois Curfman McInnes } 504c96a6f78SLois Curfman McInnes 505c96a6f78SLois Curfman McInnes #undef __FUNC__ 5065eea60f9SBarry Smith #define __FUNC__ "SNESGetSLES" /* ADIC Ignore */ 5079b94acceSBarry Smith /*@C 5089b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 5099b94acceSBarry Smith 5109b94acceSBarry Smith Input Parameter: 5119b94acceSBarry Smith . snes - the SNES context 5129b94acceSBarry Smith 5139b94acceSBarry Smith Output Parameter: 5149b94acceSBarry Smith . sles - the SLES context 5159b94acceSBarry Smith 5169b94acceSBarry Smith Notes: 5179b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 5189b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5199b94acceSBarry Smith KSP and PC contexts as well. 5209b94acceSBarry Smith 5219b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 5229b94acceSBarry Smith 5239b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 5249b94acceSBarry Smith @*/ 5259b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 5269b94acceSBarry Smith { 52777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5289b94acceSBarry Smith *sles = snes->sles; 5299b94acceSBarry Smith return 0; 5309b94acceSBarry Smith } 5319b94acceSBarry Smith /* -----------------------------------------------------------*/ 5325615d1e5SSatish Balay #undef __FUNC__ 5335615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 5349b94acceSBarry Smith /*@C 5359b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5369b94acceSBarry Smith 5379b94acceSBarry Smith Input Parameter: 5389b94acceSBarry Smith . comm - MPI communicator 5399b94acceSBarry Smith . type - type of method, one of 5409b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS 5419b94acceSBarry Smith $ (for systems of nonlinear equations) 5429b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION 5439b94acceSBarry Smith $ (for unconstrained minimization) 5449b94acceSBarry Smith 5459b94acceSBarry Smith Output Parameter: 5469b94acceSBarry Smith . outsnes - the new SNES context 5479b94acceSBarry Smith 548c1f60f51SBarry Smith Options Database Key: 54919bd6747SLois Curfman McInnes $ -snes_mf - use default matrix-free Jacobian-vector products, 55019bd6747SLois Curfman McInnes $ and no preconditioning matrix 55119bd6747SLois Curfman McInnes $ -snes_mf_operator - use default matrix-free Jacobian-vector 55219bd6747SLois Curfman McInnes $ products, and a user-provided preconditioning matrix 55319bd6747SLois Curfman McInnes $ as set by SNESSetJacobian() 55419bd6747SLois Curfman McInnes $ -snes_fd - use (slow!) finite differences to compute Jacobian 555c1f60f51SBarry Smith 5569b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5579b94acceSBarry Smith 55863a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 5599b94acceSBarry Smith @*/ 5604b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 5619b94acceSBarry Smith { 5629b94acceSBarry Smith int ierr; 5639b94acceSBarry Smith SNES snes; 5649b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 56537fcc0dbSBarry Smith 5669b94acceSBarry Smith *outsnes = 0; 5672a0acf01SLois Curfman McInnes if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS) 568d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 569f09e8eb9SSatish Balay PetscHeaderCreate(snes,_p_SNES,SNES_COOKIE,SNES_UNKNOWN_METHOD,comm); 5709b94acceSBarry Smith PLogObjectCreate(snes); 5719b94acceSBarry Smith snes->max_its = 50; 5729750a799SBarry Smith snes->max_funcs = 10000; 5739b94acceSBarry Smith snes->norm = 0.0; 5745a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 575b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 576b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5779b94acceSBarry Smith snes->atol = 1.e-10; 5785a2d0531SBarry Smith } 579b4874afaSBarry Smith else { 580b4874afaSBarry Smith snes->rtol = 1.e-8; 581b4874afaSBarry Smith snes->ttol = 0.0; 582b4874afaSBarry Smith snes->atol = 1.e-50; 583b4874afaSBarry Smith } 5849b94acceSBarry Smith snes->xtol = 1.e-8; 585b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5869b94acceSBarry Smith snes->nfuncs = 0; 5879b94acceSBarry Smith snes->nfailures = 0; 5887a00f4a9SLois Curfman McInnes snes->linear_its = 0; 589639f9d9dSBarry Smith snes->numbermonitors = 0; 5909b94acceSBarry Smith snes->data = 0; 5919b94acceSBarry Smith snes->view = 0; 5929b94acceSBarry Smith snes->computeumfunction = 0; 5939b94acceSBarry Smith snes->umfunP = 0; 5949b94acceSBarry Smith snes->fc = 0; 5959b94acceSBarry Smith snes->deltatol = 1.e-12; 5969b94acceSBarry Smith snes->fmin = -1.e30; 5979b94acceSBarry Smith snes->method_class = type; 5989b94acceSBarry Smith snes->set_method_called = 0; 599a703fe33SLois Curfman McInnes snes->setup_called = 0; 6009b94acceSBarry Smith snes->ksp_ewconv = 0; 6017d1a2b2bSBarry Smith snes->type = -1; 6026f24a144SLois Curfman McInnes snes->xmonitor = 0; 6036f24a144SLois Curfman McInnes snes->vwork = 0; 6046f24a144SLois Curfman McInnes snes->nwork = 0; 605c9005455SLois Curfman McInnes snes->conv_hist_size = 0; 606c9005455SLois Curfman McInnes snes->conv_act_size = 0; 607c9005455SLois Curfman McInnes snes->conv_hist = 0; 6089b94acceSBarry Smith 6099b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 6100452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 611eed86810SBarry Smith PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 6129b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6139b94acceSBarry Smith kctx->version = 2; 6149b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6159b94acceSBarry Smith this was too large for some test cases */ 6169b94acceSBarry Smith kctx->rtol_last = 0; 6179b94acceSBarry Smith kctx->rtol_max = .9; 6189b94acceSBarry Smith kctx->gamma = 1.0; 6199b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6209b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6219b94acceSBarry Smith kctx->threshold = .1; 6229b94acceSBarry Smith kctx->lresid_last = 0; 6239b94acceSBarry Smith kctx->norm_last = 0; 6249b94acceSBarry Smith 6259b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 6269b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 6275334005bSBarry Smith 6289b94acceSBarry Smith *outsnes = snes; 6299b94acceSBarry Smith return 0; 6309b94acceSBarry Smith } 6319b94acceSBarry Smith 6329b94acceSBarry Smith /* --------------------------------------------------------------- */ 6335615d1e5SSatish Balay #undef __FUNC__ 6345eea60f9SBarry Smith #define __FUNC__ "SNESSetFunction" /* ADIC Ignore */ 6359b94acceSBarry Smith /*@C 6369b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6379b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6389b94acceSBarry Smith equations. 6399b94acceSBarry Smith 6409b94acceSBarry Smith Input Parameters: 6419b94acceSBarry Smith . snes - the SNES context 6429b94acceSBarry Smith . func - function evaluation routine 6439b94acceSBarry Smith . r - vector to store function value 6442cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 6452cd2dfdcSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6469b94acceSBarry Smith 6479b94acceSBarry Smith Calling sequence of func: 648313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Vec f,void *ctx); 6499b94acceSBarry Smith 6509b94acceSBarry Smith . x - input vector 651313e4042SLois Curfman McInnes . f - function vector 6522cd2dfdcSLois Curfman McInnes . ctx - optional user-defined function context 6539b94acceSBarry Smith 6549b94acceSBarry Smith Notes: 6559b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6569b94acceSBarry Smith $ f'(x) x = -f(x), 6579b94acceSBarry Smith $ where f'(x) denotes the Jacobian matrix and f(x) is the function. 6589b94acceSBarry Smith 6599b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6609b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6619b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 6629b94acceSBarry Smith 6639b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6649b94acceSBarry Smith 665a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6669b94acceSBarry Smith @*/ 6675334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 6689b94acceSBarry Smith { 66977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 670e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 6719b94acceSBarry Smith snes->computefunction = func; 6729b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6739b94acceSBarry Smith snes->funP = ctx; 6749b94acceSBarry Smith return 0; 6759b94acceSBarry Smith } 6769b94acceSBarry Smith 6775615d1e5SSatish Balay #undef __FUNC__ 6785615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 6799b94acceSBarry Smith /*@ 6809b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6819b94acceSBarry Smith SNESSetFunction(). 6829b94acceSBarry Smith 6839b94acceSBarry Smith Input Parameters: 6849b94acceSBarry Smith . snes - the SNES context 6859b94acceSBarry Smith . x - input vector 6869b94acceSBarry Smith 6879b94acceSBarry Smith Output Parameter: 6883638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6899b94acceSBarry Smith 6901bffabb2SLois Curfman McInnes Notes: 6919b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6929b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6939b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6949b94acceSBarry Smith 6959b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6969b94acceSBarry Smith 697a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6989b94acceSBarry Smith @*/ 6999b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 7009b94acceSBarry Smith { 7019b94acceSBarry Smith int ierr; 7029b94acceSBarry Smith 70374679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 704e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 7059b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 7069b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 707ae3c334cSLois Curfman McInnes snes->nfuncs++; 7089b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 7099b94acceSBarry Smith return 0; 7109b94acceSBarry Smith } 7119b94acceSBarry Smith 7125615d1e5SSatish Balay #undef __FUNC__ 7135eea60f9SBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" /* ADIC Ignore */ 7149b94acceSBarry Smith /*@C 7159b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 7169b94acceSBarry Smith unconstrained minimization. 7179b94acceSBarry Smith 7189b94acceSBarry Smith Input Parameters: 7199b94acceSBarry Smith . snes - the SNES context 7209b94acceSBarry Smith . func - function evaluation routine 721044dda88SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 722044dda88SLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 7239b94acceSBarry Smith 7249b94acceSBarry Smith Calling sequence of func: 7259b94acceSBarry Smith . func (SNES snes,Vec x,double *f,void *ctx); 7269b94acceSBarry Smith 7279b94acceSBarry Smith . x - input vector 7289b94acceSBarry Smith . f - function 729044dda88SLois Curfman McInnes . ctx - [optional] user-defined function context 7309b94acceSBarry Smith 7319b94acceSBarry Smith Notes: 7329b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7339b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7349b94acceSBarry Smith SNESSetFunction(). 7359b94acceSBarry Smith 7369b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 7379b94acceSBarry Smith 738a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 739a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 7409b94acceSBarry Smith @*/ 7419b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 7429b94acceSBarry Smith void *ctx) 7439b94acceSBarry Smith { 74477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 745e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7469b94acceSBarry Smith snes->computeumfunction = func; 7479b94acceSBarry Smith snes->umfunP = ctx; 7489b94acceSBarry Smith return 0; 7499b94acceSBarry Smith } 7509b94acceSBarry Smith 7515615d1e5SSatish Balay #undef __FUNC__ 7525615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 7539b94acceSBarry Smith /*@ 7549b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 7559b94acceSBarry Smith set with SNESSetMinimizationFunction(). 7569b94acceSBarry Smith 7579b94acceSBarry Smith Input Parameters: 7589b94acceSBarry Smith . snes - the SNES context 7599b94acceSBarry Smith . x - input vector 7609b94acceSBarry Smith 7619b94acceSBarry Smith Output Parameter: 7629b94acceSBarry Smith . y - function value 7639b94acceSBarry Smith 7649b94acceSBarry Smith Notes: 7659b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 7669b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7679b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 768a86d99e1SLois Curfman McInnes 769a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 770a86d99e1SLois Curfman McInnes 771a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 772a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 7739b94acceSBarry Smith @*/ 7749b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 7759b94acceSBarry Smith { 7769b94acceSBarry Smith int ierr; 777e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 7789b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 7799b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 780ae3c334cSLois Curfman McInnes snes->nfuncs++; 7819b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7829b94acceSBarry Smith return 0; 7839b94acceSBarry Smith } 7849b94acceSBarry Smith 7855615d1e5SSatish Balay #undef __FUNC__ 7865eea60f9SBarry Smith #define __FUNC__ "SNESSetGradient" /* ADIC Ignore */ 7879b94acceSBarry Smith /*@C 7889b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 7899b94acceSBarry Smith vector for use by the SNES routines. 7909b94acceSBarry Smith 7919b94acceSBarry Smith Input Parameters: 7929b94acceSBarry Smith . snes - the SNES context 7939b94acceSBarry Smith . func - function evaluation routine 794044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 795044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 7969b94acceSBarry Smith . r - vector to store gradient value 7979b94acceSBarry Smith 7989b94acceSBarry Smith Calling sequence of func: 7999b94acceSBarry Smith . func (SNES, Vec x, Vec g, void *ctx); 8009b94acceSBarry Smith 8019b94acceSBarry Smith . x - input vector 8029b94acceSBarry Smith . g - gradient vector 803044dda88SLois Curfman McInnes . ctx - optional user-defined gradient context 8049b94acceSBarry Smith 8059b94acceSBarry Smith Notes: 8069b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 8079b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 8089b94acceSBarry Smith SNESSetFunction(). 8099b94acceSBarry Smith 8109b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 8119b94acceSBarry Smith 812a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 813a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 8149b94acceSBarry Smith @*/ 81574679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 8169b94acceSBarry Smith { 81777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 818e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8199b94acceSBarry Smith snes->computefunction = func; 8209b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 8219b94acceSBarry Smith snes->funP = ctx; 8229b94acceSBarry Smith return 0; 8239b94acceSBarry Smith } 8249b94acceSBarry Smith 8255615d1e5SSatish Balay #undef __FUNC__ 8265615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 8279b94acceSBarry Smith /*@ 828a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 829a86d99e1SLois Curfman McInnes SNESSetGradient(). 8309b94acceSBarry Smith 8319b94acceSBarry Smith Input Parameters: 8329b94acceSBarry Smith . snes - the SNES context 8339b94acceSBarry Smith . x - input vector 8349b94acceSBarry Smith 8359b94acceSBarry Smith Output Parameter: 8369b94acceSBarry Smith . y - gradient vector 8379b94acceSBarry Smith 8389b94acceSBarry Smith Notes: 8399b94acceSBarry Smith SNESComputeGradient() is valid only for 8409b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 8419b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 842a86d99e1SLois Curfman McInnes 843a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 844a86d99e1SLois Curfman McInnes 845a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 846a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 8479b94acceSBarry Smith @*/ 8489b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 8499b94acceSBarry Smith { 8509b94acceSBarry Smith int ierr; 85174679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 852e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8539b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 8549b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 8559b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 8569b94acceSBarry Smith return 0; 8579b94acceSBarry Smith } 8589b94acceSBarry Smith 8595615d1e5SSatish Balay #undef __FUNC__ 8605615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 86162fef451SLois Curfman McInnes /*@ 86262fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 86362fef451SLois Curfman McInnes set with SNESSetJacobian(). 86462fef451SLois Curfman McInnes 86562fef451SLois Curfman McInnes Input Parameters: 86662fef451SLois Curfman McInnes . snes - the SNES context 86762fef451SLois Curfman McInnes . x - input vector 86862fef451SLois Curfman McInnes 86962fef451SLois Curfman McInnes Output Parameters: 87062fef451SLois Curfman McInnes . A - Jacobian matrix 87162fef451SLois Curfman McInnes . B - optional preconditioning matrix 87262fef451SLois Curfman McInnes . flag - flag indicating matrix structure 87362fef451SLois Curfman McInnes 87462fef451SLois Curfman McInnes Notes: 87562fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 87662fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 87762fef451SLois Curfman McInnes 878dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 879dc5a77f8SLois Curfman McInnes flag parameter. 88062fef451SLois Curfman McInnes 88162fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 88262fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 88362fef451SLois Curfman McInnes methods is SNESComputeJacobian(). 88462fef451SLois Curfman McInnes 88562fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 88662fef451SLois Curfman McInnes 88762fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 88862fef451SLois Curfman McInnes @*/ 8899b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8909b94acceSBarry Smith { 8919b94acceSBarry Smith int ierr; 89274679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 893e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 8949b94acceSBarry Smith if (!snes->computejacobian) return 0; 8959b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 896c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 8979b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 8989b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 8996d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 90077c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 90177c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9029b94acceSBarry Smith return 0; 9039b94acceSBarry Smith } 9049b94acceSBarry Smith 9055615d1e5SSatish Balay #undef __FUNC__ 9065615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 90762fef451SLois Curfman McInnes /*@ 90862fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 90962fef451SLois Curfman McInnes set with SNESSetHessian(). 91062fef451SLois Curfman McInnes 91162fef451SLois Curfman McInnes Input Parameters: 91262fef451SLois Curfman McInnes . snes - the SNES context 91362fef451SLois Curfman McInnes . x - input vector 91462fef451SLois Curfman McInnes 91562fef451SLois Curfman McInnes Output Parameters: 91662fef451SLois Curfman McInnes . A - Hessian matrix 91762fef451SLois Curfman McInnes . B - optional preconditioning matrix 91862fef451SLois Curfman McInnes . flag - flag indicating matrix structure 91962fef451SLois Curfman McInnes 92062fef451SLois Curfman McInnes Notes: 92162fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 92262fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 92362fef451SLois Curfman McInnes 924dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 925dc5a77f8SLois Curfman McInnes flag parameter. 92662fef451SLois Curfman McInnes 92762fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 92862fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 92962fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 93062fef451SLois Curfman McInnes 93162fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 93262fef451SLois Curfman McInnes 933a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 934a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 93562fef451SLois Curfman McInnes @*/ 93662fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 9379b94acceSBarry Smith { 9389b94acceSBarry Smith int ierr; 93974679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 940e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 9419b94acceSBarry Smith if (!snes->computejacobian) return 0; 94262fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 9430f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 94462fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 94562fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 9460f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 94777c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 94877c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9499b94acceSBarry Smith return 0; 9509b94acceSBarry Smith } 9519b94acceSBarry Smith 9525615d1e5SSatish Balay #undef __FUNC__ 9535eea60f9SBarry Smith #define __FUNC__ "SNESSetJacobian" /* ADIC Ignore */ 9549b94acceSBarry Smith /*@C 9559b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 956044dda88SLois Curfman McInnes location to store the matrix. 9579b94acceSBarry Smith 9589b94acceSBarry Smith Input Parameters: 9599b94acceSBarry Smith . snes - the SNES context 9609b94acceSBarry Smith . A - Jacobian matrix 9619b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 9629b94acceSBarry Smith . func - Jacobian evaluation routine 9632cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 9642cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 9659b94acceSBarry Smith 9669b94acceSBarry Smith Calling sequence of func: 967313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 9689b94acceSBarry Smith 9699b94acceSBarry Smith . x - input vector 9709b94acceSBarry Smith . A - Jacobian matrix 9719b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 972ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 973ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 9742cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Jacobian context 9759b94acceSBarry Smith 9769b94acceSBarry Smith Notes: 977dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 9782cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 979ac21db08SLois Curfman McInnes 980ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 9819b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 9829b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 9839b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 9849b94acceSBarry Smith throughout the global iterations. 9859b94acceSBarry Smith 9869b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 9879b94acceSBarry Smith 988ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 9899b94acceSBarry Smith @*/ 9909b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 9919b94acceSBarry Smith MatStructure*,void*),void *ctx) 9929b94acceSBarry Smith { 99377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 99474679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 995e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 9969b94acceSBarry Smith snes->computejacobian = func; 9979b94acceSBarry Smith snes->jacP = ctx; 9989b94acceSBarry Smith snes->jacobian = A; 9999b94acceSBarry Smith snes->jacobian_pre = B; 10009b94acceSBarry Smith return 0; 10019b94acceSBarry Smith } 100262fef451SLois Curfman McInnes 10035615d1e5SSatish Balay #undef __FUNC__ 10045eea60f9SBarry Smith #define __FUNC__ "SNESGetJacobian" /* ADIC Ignore */ 1005b4fd4287SBarry Smith /*@ 1006b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1007b4fd4287SBarry Smith provided context for evaluating the Jacobian. 1008b4fd4287SBarry Smith 1009b4fd4287SBarry Smith Input Parameter: 1010b4fd4287SBarry Smith . snes - the nonlinear solver context 1011b4fd4287SBarry Smith 1012b4fd4287SBarry Smith Output Parameters: 1013b4fd4287SBarry Smith . A - location to stash Jacobian matrix (or PETSC_NULL) 1014b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 1015b4fd4287SBarry Smith . ctx - location to stash Jacobian ctx (or PETSC_NULL) 1016b4fd4287SBarry Smith 1017b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 1018b4fd4287SBarry Smith @*/ 1019b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1020b4fd4287SBarry Smith { 102174679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) 1022e3372554SBarry Smith SETERRQ(1,0,"For SNES_NONLINEAR_EQUATIONS only"); 1023b4fd4287SBarry Smith if (A) *A = snes->jacobian; 1024b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 1025b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 1026b4fd4287SBarry Smith return 0; 1027b4fd4287SBarry Smith } 1028b4fd4287SBarry Smith 10295615d1e5SSatish Balay #undef __FUNC__ 10305eea60f9SBarry Smith #define __FUNC__ "SNESSetHessian" /* ADIC Ignore */ 10319b94acceSBarry Smith /*@C 10329b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 1033044dda88SLois Curfman McInnes location to store the matrix. 10349b94acceSBarry Smith 10359b94acceSBarry Smith Input Parameters: 10369b94acceSBarry Smith . snes - the SNES context 10379b94acceSBarry Smith . A - Hessian matrix 10389b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 10399b94acceSBarry Smith . func - Jacobian evaluation routine 1040313e4042SLois Curfman McInnes . ctx - [optional] user-defined context for private data for the 1041313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 10429b94acceSBarry Smith 10439b94acceSBarry Smith Calling sequence of func: 1044313e4042SLois Curfman McInnes . func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 10459b94acceSBarry Smith 10469b94acceSBarry Smith . x - input vector 10479b94acceSBarry Smith . A - Hessian matrix 10489b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1049ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1050ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 10512cd2dfdcSLois Curfman McInnes . ctx - [optional] user-defined Hessian context 10529b94acceSBarry Smith 10539b94acceSBarry Smith Notes: 1054dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10552cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1056ac21db08SLois Curfman McInnes 10579b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 10589b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 10599b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10609b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10619b94acceSBarry Smith throughout the global iterations. 10629b94acceSBarry Smith 10639b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 10649b94acceSBarry Smith 1065ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 10669b94acceSBarry Smith @*/ 10679b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 10689b94acceSBarry Smith MatStructure*,void*),void *ctx) 10699b94acceSBarry Smith { 107077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 107174679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1072e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 10739b94acceSBarry Smith snes->computejacobian = func; 10749b94acceSBarry Smith snes->jacP = ctx; 10759b94acceSBarry Smith snes->jacobian = A; 10769b94acceSBarry Smith snes->jacobian_pre = B; 10779b94acceSBarry Smith return 0; 10789b94acceSBarry Smith } 10799b94acceSBarry Smith 10805615d1e5SSatish Balay #undef __FUNC__ 10815eea60f9SBarry Smith #define __FUNC__ "SNESGetHessian" /* ADIC Ignore */ 108262fef451SLois Curfman McInnes /*@ 108362fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 108462fef451SLois Curfman McInnes provided context for evaluating the Hessian. 108562fef451SLois Curfman McInnes 108662fef451SLois Curfman McInnes Input Parameter: 108762fef451SLois Curfman McInnes . snes - the nonlinear solver context 108862fef451SLois Curfman McInnes 108962fef451SLois Curfman McInnes Output Parameters: 109062fef451SLois Curfman McInnes . A - location to stash Hessian matrix (or PETSC_NULL) 109162fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 109262fef451SLois Curfman McInnes . ctx - location to stash Hessian ctx (or PETSC_NULL) 109362fef451SLois Curfman McInnes 109462fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 109562fef451SLois Curfman McInnes @*/ 109662fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 109762fef451SLois Curfman McInnes { 109874679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) 1099e3372554SBarry Smith SETERRQ(1,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 110062fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 110162fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 110262fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 110362fef451SLois Curfman McInnes return 0; 110462fef451SLois Curfman McInnes } 110562fef451SLois Curfman McInnes 11069b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 11079b94acceSBarry Smith 11085615d1e5SSatish Balay #undef __FUNC__ 11095615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 11109b94acceSBarry Smith /*@ 11119b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1112272ac6f2SLois Curfman McInnes of a nonlinear solver. 11139b94acceSBarry Smith 11149b94acceSBarry Smith Input Parameter: 11159b94acceSBarry Smith . snes - the SNES context 11168ddd3da0SLois Curfman McInnes . x - the solution vector 11179b94acceSBarry Smith 1118272ac6f2SLois Curfman McInnes Notes: 1119272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1120272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1121272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1122272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1123272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1124272ac6f2SLois Curfman McInnes 11259b94acceSBarry Smith .keywords: SNES, nonlinear, setup 11269b94acceSBarry Smith 11279b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 11289b94acceSBarry Smith @*/ 11298ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 11309b94acceSBarry Smith { 1131272ac6f2SLois Curfman McInnes int ierr, flg; 113277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 113377c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 11348ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 11359b94acceSBarry Smith 1136c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1137c1f60f51SBarry Smith /* 1138c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1139dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1140c1f60f51SBarry Smith */ 1141c1f60f51SBarry Smith if (flg) { 1142c1f60f51SBarry Smith Mat J; 1143c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1144c1f60f51SBarry Smith PLogObjectParent(snes,J); 1145c1f60f51SBarry Smith snes->mfshell = J; 1146c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1147c1f60f51SBarry Smith snes->jacobian = J; 1148c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1149c1f60f51SBarry Smith } 1150c1f60f51SBarry Smith else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1151c1f60f51SBarry Smith snes->jacobian = J; 1152c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1153e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free operator option"); 1154c1f60f51SBarry Smith } 1155272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1156c1f60f51SBarry Smith /* 1157dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1158c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1159c1f60f51SBarry Smith */ 116031615d04SLois Curfman McInnes if (flg) { 1161272ac6f2SLois Curfman McInnes Mat J; 1162272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1163272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1164272ac6f2SLois Curfman McInnes snes->mfshell = J; 116583e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 116683e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 116794a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 116883e56afdSLois Curfman McInnes } 116983e56afdSLois Curfman McInnes else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 117083e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 117194a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1172e3372554SBarry Smith } else SETERRQ(1,0,"Method class doesn't support matrix-free option"); 1173272ac6f2SLois Curfman McInnes } 11749b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1175e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1176e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetFunction() first"); 1177e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetJacobian() first"); 1178e3372554SBarry Smith if (snes->vec_func == snes->vec_sol) SETERRQ(1,0,"Solution vector cannot be function vector"); 1179a703fe33SLois Curfman McInnes 1180a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 118140191667SLois Curfman McInnes if (snes->ksp_ewconv && snes->type != SNES_EQ_TR) { 1182a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1183a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1184a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1185a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1186a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1187a703fe33SLois Curfman McInnes } 11889b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1189e3372554SBarry Smith if (!snes->vec_func) SETERRQ(1,0,"Must call SNESSetGradient() first"); 1190e3372554SBarry Smith if (!snes->computefunction) SETERRQ(1,0,"Must call SNESSetGradient() first"); 119137fcc0dbSBarry Smith if (!snes->computeumfunction) 1192e3372554SBarry Smith SETERRQ(1,0,"Must call SNESSetMinimizationFunction() first"); 1193e3372554SBarry Smith if (!snes->jacobian) SETERRQ(1,0,"Must call SNESSetHessian() first"); 1194e3372554SBarry Smith } else SETERRQ(1,0,"Unknown method class"); 1195a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 1196a703fe33SLois Curfman McInnes snes->setup_called = 1; 1197a703fe33SLois Curfman McInnes return 0; 11989b94acceSBarry Smith } 11999b94acceSBarry Smith 12005615d1e5SSatish Balay #undef __FUNC__ 12015eea60f9SBarry Smith #define __FUNC__ "SNESDestroy" /* ADIC Ignore */ 12029b94acceSBarry Smith /*@C 12039b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 12049b94acceSBarry Smith with SNESCreate(). 12059b94acceSBarry Smith 12069b94acceSBarry Smith Input Parameter: 12079b94acceSBarry Smith . snes - the SNES context 12089b94acceSBarry Smith 12099b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 12109b94acceSBarry Smith 121163a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 12129b94acceSBarry Smith @*/ 12139b94acceSBarry Smith int SNESDestroy(SNES snes) 12149b94acceSBarry Smith { 12159b94acceSBarry Smith int ierr; 121677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12179750a799SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)((PetscObject)snes); CHKERRQ(ierr);} 12180452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 12199b94acceSBarry Smith if (snes->mfshell) MatDestroy(snes->mfshell); 12209b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 12213e2e452bSBarry Smith if (snes->xmonitor) SNESLGMonitorDestroy(snes->xmonitor); 12226f24a144SLois Curfman McInnes if (snes->vwork) VecDestroyVecs(snes->vwork,snes->nvwork); 12239b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 12240452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 12259b94acceSBarry Smith return 0; 12269b94acceSBarry Smith } 12279b94acceSBarry Smith 12289b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 12299b94acceSBarry Smith 12305615d1e5SSatish Balay #undef __FUNC__ 12315615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 12329b94acceSBarry Smith /*@ 1233d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 12349b94acceSBarry Smith 12359b94acceSBarry Smith Input Parameters: 12369b94acceSBarry Smith . snes - the SNES context 123733174efeSLois Curfman McInnes . atol - absolute convergence tolerance 123833174efeSLois Curfman McInnes . rtol - relative convergence tolerance 123933174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 124033174efeSLois Curfman McInnes of the change in the solution between steps 124133174efeSLois Curfman McInnes . maxit - maximum number of iterations 124233174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 12439b94acceSBarry Smith 124433174efeSLois Curfman McInnes Options Database Keys: 124533174efeSLois Curfman McInnes $ -snes_atol <atol> 124633174efeSLois Curfman McInnes $ -snes_rtol <rtol> 124733174efeSLois Curfman McInnes $ -snes_stol <stol> 124833174efeSLois Curfman McInnes $ -snes_max_it <maxit> 124933174efeSLois Curfman McInnes $ -snes_max_funcs <maxf> 12509b94acceSBarry Smith 1251d7a720efSLois Curfman McInnes Notes: 12529b94acceSBarry Smith The default maximum number of iterations is 50. 12539b94acceSBarry Smith The default maximum number of function evaluations is 1000. 12549b94acceSBarry Smith 125533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 12569b94acceSBarry Smith 1257d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 12589b94acceSBarry Smith @*/ 1259d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 12609b94acceSBarry Smith { 126177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1262d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1263d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1264d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1265d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1266d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 12679b94acceSBarry Smith return 0; 12689b94acceSBarry Smith } 12699b94acceSBarry Smith 12705615d1e5SSatish Balay #undef __FUNC__ 12715615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 12729b94acceSBarry Smith /*@ 127333174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 127433174efeSLois Curfman McInnes 127533174efeSLois Curfman McInnes Input Parameters: 127633174efeSLois Curfman McInnes . snes - the SNES context 127733174efeSLois Curfman McInnes . atol - absolute convergence tolerance 127833174efeSLois Curfman McInnes . rtol - relative convergence tolerance 127933174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 128033174efeSLois Curfman McInnes of the change in the solution between steps 128133174efeSLois Curfman McInnes . maxit - maximum number of iterations 128233174efeSLois Curfman McInnes . maxf - maximum number of function evaluations 128333174efeSLois Curfman McInnes 128433174efeSLois Curfman McInnes Notes: 128533174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 128633174efeSLois Curfman McInnes 128733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 128833174efeSLois Curfman McInnes 128933174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 129033174efeSLois Curfman McInnes @*/ 129133174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 129233174efeSLois Curfman McInnes { 129333174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 129433174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 129533174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 129633174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 129733174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 129833174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 129933174efeSLois Curfman McInnes return 0; 130033174efeSLois Curfman McInnes } 130133174efeSLois Curfman McInnes 13025615d1e5SSatish Balay #undef __FUNC__ 13035615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 130433174efeSLois Curfman McInnes /*@ 13059b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 13069b94acceSBarry Smith 13079b94acceSBarry Smith Input Parameters: 13089b94acceSBarry Smith . snes - the SNES context 13099b94acceSBarry Smith . tol - tolerance 13109b94acceSBarry Smith 13119b94acceSBarry Smith Options Database Key: 1312d7a720efSLois Curfman McInnes $ -snes_trtol <tol> 13139b94acceSBarry Smith 13149b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 13159b94acceSBarry Smith 1316d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 13179b94acceSBarry Smith @*/ 13189b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 13199b94acceSBarry Smith { 132077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13219b94acceSBarry Smith snes->deltatol = tol; 13229b94acceSBarry Smith return 0; 13239b94acceSBarry Smith } 13249b94acceSBarry Smith 13255615d1e5SSatish Balay #undef __FUNC__ 13265615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 13279b94acceSBarry Smith /*@ 132877c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 13299b94acceSBarry Smith for unconstrained minimization solvers. 13309b94acceSBarry Smith 13319b94acceSBarry Smith Input Parameters: 13329b94acceSBarry Smith . snes - the SNES context 13339b94acceSBarry Smith . ftol - minimum function tolerance 13349b94acceSBarry Smith 13359b94acceSBarry Smith Options Database Key: 1336d7a720efSLois Curfman McInnes $ -snes_fmin <ftol> 13379b94acceSBarry Smith 13389b94acceSBarry Smith Note: 133977c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 13409b94acceSBarry Smith methods only. 13419b94acceSBarry Smith 13429b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 13439b94acceSBarry Smith 1344d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 13459b94acceSBarry Smith @*/ 134677c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 13479b94acceSBarry Smith { 134877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13499b94acceSBarry Smith snes->fmin = ftol; 13509b94acceSBarry Smith return 0; 13519b94acceSBarry Smith } 13529b94acceSBarry Smith 13539b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 13549b94acceSBarry Smith 13555615d1e5SSatish Balay #undef __FUNC__ 13565eea60f9SBarry Smith #define __FUNC__ "SNESSetMonitor" /* ADIC Ignore */ 13579b94acceSBarry Smith /*@C 13589b94acceSBarry Smith SNESSetMonitor - Sets the function that is to be used at every 13599b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 13609b94acceSBarry Smith progress. 13619b94acceSBarry Smith 13629b94acceSBarry Smith Input Parameters: 13639b94acceSBarry Smith . snes - the SNES context 13649b94acceSBarry Smith . func - monitoring routine 1365044dda88SLois Curfman McInnes . mctx - [optional] user-defined context for private data for the 1366044dda88SLois Curfman McInnes monitor routine (may be PETSC_NULL) 13679b94acceSBarry Smith 13689b94acceSBarry Smith Calling sequence of func: 1369682d7d0cSBarry Smith int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 13709b94acceSBarry Smith 13719b94acceSBarry Smith $ snes - the SNES context 13729b94acceSBarry Smith $ its - iteration number 1373044dda88SLois Curfman McInnes $ mctx - [optional] monitoring context 13749b94acceSBarry Smith $ 13759b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 13769b94acceSBarry Smith $ norm - 2-norm function value (may be estimated) 13779b94acceSBarry Smith $ 13789b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 13799b94acceSBarry Smith $ norm - 2-norm gradient value (may be estimated) 13809b94acceSBarry Smith 13819665c990SLois Curfman McInnes Options Database Keys: 13829665c990SLois Curfman McInnes $ -snes_monitor : sets SNESDefaultMonitor() 13839665c990SLois Curfman McInnes $ -snes_xmonitor : sets line graph monitor, 13849665c990SLois Curfman McInnes $ uses SNESLGMonitorCreate() 13859665c990SLois Curfman McInnes $ -snes_cancelmonitors : cancels all monitors that have 13869665c990SLois Curfman McInnes $ been hardwired into a code by 13879665c990SLois Curfman McInnes $ calls to SNESSetMonitor(), but 13889665c990SLois Curfman McInnes $ does not cancel those set via 13899665c990SLois Curfman McInnes $ the options database. 13909665c990SLois Curfman McInnes 13919665c990SLois Curfman McInnes 1392639f9d9dSBarry Smith Notes: 13936bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 13946bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 13956bc08f3fSLois Curfman McInnes order in which they were set. 1396639f9d9dSBarry Smith 13979b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 13989b94acceSBarry Smith 13999b94acceSBarry Smith .seealso: SNESDefaultMonitor() 14009b94acceSBarry Smith @*/ 140174679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 14029b94acceSBarry Smith { 1403639f9d9dSBarry Smith if (!func) { 1404639f9d9dSBarry Smith snes->numbermonitors = 0; 1405639f9d9dSBarry Smith return 0; 1406639f9d9dSBarry Smith } 1407639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1408e3372554SBarry Smith SETERRQ(1,0,"Too many monitors set"); 1409639f9d9dSBarry Smith } 1410639f9d9dSBarry Smith 1411639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1412639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 14139b94acceSBarry Smith return 0; 14149b94acceSBarry Smith } 14159b94acceSBarry Smith 14165615d1e5SSatish Balay #undef __FUNC__ 14175eea60f9SBarry Smith #define __FUNC__ "SNESSetConvergenceTest" /* ADIC Ignore */ 14189b94acceSBarry Smith /*@C 14199b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 14209b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 14219b94acceSBarry Smith 14229b94acceSBarry Smith Input Parameters: 14239b94acceSBarry Smith . snes - the SNES context 14249b94acceSBarry Smith . func - routine to test for convergence 1425044dda88SLois Curfman McInnes . cctx - [optional] context for private data for the convergence routine 1426044dda88SLois Curfman McInnes (may be PETSC_NULL) 14279b94acceSBarry Smith 14289b94acceSBarry Smith Calling sequence of func: 14299b94acceSBarry Smith int func (SNES snes,double xnorm,double gnorm, 14309b94acceSBarry Smith double f,void *cctx) 14319b94acceSBarry Smith 14329b94acceSBarry Smith $ snes - the SNES context 1433044dda88SLois Curfman McInnes $ cctx - [optional] convergence context 14349b94acceSBarry Smith $ xnorm - 2-norm of current iterate 14359b94acceSBarry Smith $ 14369b94acceSBarry Smith $ SNES_NONLINEAR_EQUATIONS methods: 14379b94acceSBarry Smith $ gnorm - 2-norm of current step 14389b94acceSBarry Smith $ f - 2-norm of function 14399b94acceSBarry Smith $ 14409b94acceSBarry Smith $ SNES_UNCONSTRAINED_MINIMIZATION methods: 14419b94acceSBarry Smith $ gnorm - 2-norm of current gradient 14429b94acceSBarry Smith $ f - function value 14439b94acceSBarry Smith 14449b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 14459b94acceSBarry Smith 144640191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 144740191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 14489b94acceSBarry Smith @*/ 144974679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 14509b94acceSBarry Smith { 14519b94acceSBarry Smith (snes)->converged = func; 14529b94acceSBarry Smith (snes)->cnvP = cctx; 14539b94acceSBarry Smith return 0; 14549b94acceSBarry Smith } 14559b94acceSBarry Smith 14565615d1e5SSatish Balay #undef __FUNC__ 14575eea60f9SBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" /* ADIC Ignore */ 1458c9005455SLois Curfman McInnes /*@ 1459c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1460c9005455SLois Curfman McInnes 1461c9005455SLois Curfman McInnes Input Parameters: 1462c9005455SLois Curfman McInnes . snes - iterative context obtained from SNESCreate() 1463c9005455SLois Curfman McInnes . a - array to hold history 1464c9005455SLois Curfman McInnes . na - size of a 1465c9005455SLois Curfman McInnes 1466c9005455SLois Curfman McInnes Notes: 1467c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1468c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1469c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1470c9005455SLois Curfman McInnes at each step. 1471c9005455SLois Curfman McInnes 1472c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1473c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1474c9005455SLois Curfman McInnes during the section of code that is being timed. 1475c9005455SLois Curfman McInnes 1476c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1477c9005455SLois Curfman McInnes @*/ 1478c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na) 1479c9005455SLois Curfman McInnes { 1480c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1481c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1482c9005455SLois Curfman McInnes snes->conv_hist = a; 1483c9005455SLois Curfman McInnes snes->conv_hist_size = na; 1484c9005455SLois Curfman McInnes return 0; 1485c9005455SLois Curfman McInnes } 1486c9005455SLois Curfman McInnes 1487c9005455SLois Curfman McInnes #undef __FUNC__ 14885615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 14899b94acceSBarry Smith /* 14909b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 14919b94acceSBarry Smith positive parameter delta. 14929b94acceSBarry Smith 14939b94acceSBarry Smith Input Parameters: 14949b94acceSBarry Smith . snes - the SNES context 14959b94acceSBarry Smith . y - approximate solution of linear system 14969b94acceSBarry Smith . fnorm - 2-norm of current function 14979b94acceSBarry Smith . delta - trust region size 14989b94acceSBarry Smith 14999b94acceSBarry Smith Output Parameters: 15009b94acceSBarry Smith . gpnorm - predicted function norm at the new point, assuming local 15019b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 15029b94acceSBarry Smith region, and exceeds zero otherwise. 15039b94acceSBarry Smith . ynorm - 2-norm of the step 15049b94acceSBarry Smith 15059b94acceSBarry Smith Note: 150640191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 15079b94acceSBarry Smith is set to be the maximum allowable step size. 15089b94acceSBarry Smith 15099b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 15109b94acceSBarry Smith */ 15119b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 15129b94acceSBarry Smith double *gpnorm,double *ynorm) 15139b94acceSBarry Smith { 15149b94acceSBarry Smith double norm; 15159b94acceSBarry Smith Scalar cnorm; 1516cddf8d76SBarry Smith VecNorm(y,NORM_2, &norm ); 15179b94acceSBarry Smith if (norm > *delta) { 15189b94acceSBarry Smith norm = *delta/norm; 15199b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 15209b94acceSBarry Smith cnorm = norm; 15219b94acceSBarry Smith VecScale( &cnorm, y ); 15229b94acceSBarry Smith *ynorm = *delta; 15239b94acceSBarry Smith } else { 15249b94acceSBarry Smith *gpnorm = 0.0; 15259b94acceSBarry Smith *ynorm = norm; 15269b94acceSBarry Smith } 15279b94acceSBarry Smith return 0; 15289b94acceSBarry Smith } 15299b94acceSBarry Smith 15305615d1e5SSatish Balay #undef __FUNC__ 15315615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 15329b94acceSBarry Smith /*@ 15339b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 153463a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 15359b94acceSBarry Smith 15369b94acceSBarry Smith Input Parameter: 15379b94acceSBarry Smith . snes - the SNES context 15388ddd3da0SLois Curfman McInnes . x - the solution vector 15399b94acceSBarry Smith 15409b94acceSBarry Smith Output Parameter: 15419b94acceSBarry Smith its - number of iterations until termination 15429b94acceSBarry Smith 15438ddd3da0SLois Curfman McInnes Note: 15448ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 15458ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 15468ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 15478ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 15488ddd3da0SLois Curfman McInnes 15499b94acceSBarry Smith .keywords: SNES, nonlinear, solve 15509b94acceSBarry Smith 155163a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 15529b94acceSBarry Smith @*/ 15538ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 15549b94acceSBarry Smith { 15553c7409f5SSatish Balay int ierr, flg; 1556052efed2SBarry Smith 155777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 155874679c65SBarry Smith PetscValidIntPointer(its); 1559c4fc05e7SBarry Smith if (!snes->setup_called) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1560c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 15619b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 1562c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 15639b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 15649b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 15653c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 15666d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 15679b94acceSBarry Smith return 0; 15689b94acceSBarry Smith } 15699b94acceSBarry Smith 15709b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 157137fcc0dbSBarry Smith static NRList *__SNESList = 0; 15729b94acceSBarry Smith 15735615d1e5SSatish Balay #undef __FUNC__ 15745615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 15759b94acceSBarry Smith /*@ 15764b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 15779b94acceSBarry Smith 15789b94acceSBarry Smith Input Parameters: 15799b94acceSBarry Smith . snes - the SNES context 15809b94acceSBarry Smith . method - a known method 15819b94acceSBarry Smith 1582ae12b187SLois Curfman McInnes Options Database Command: 1583ae12b187SLois Curfman McInnes $ -snes_type <method> 1584ae12b187SLois Curfman McInnes $ Use -help for a list of available methods 1585ae12b187SLois Curfman McInnes $ (for instance, ls or tr) 1586ae12b187SLois Curfman McInnes 15879b94acceSBarry Smith Notes: 15889b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 15899b94acceSBarry Smith $ Systems of nonlinear equations: 159040191667SLois Curfman McInnes $ SNES_EQ_LS - Newton's method with line search 159140191667SLois Curfman McInnes $ SNES_EQ_TR - Newton's method with trust region 15929b94acceSBarry Smith $ Unconstrained minimization: 159340191667SLois Curfman McInnes $ SNES_UM_TR - Newton's method with trust region 159440191667SLois Curfman McInnes $ SNES_UM_LS - Newton's method with line search 15959b94acceSBarry Smith 1596ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1597ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1598ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1599ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1600ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1601ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1602ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1603ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1604ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1605ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1606a703fe33SLois Curfman McInnes 1607f536c699SSatish Balay .keywords: SNES, set, method 16089b94acceSBarry Smith @*/ 16094b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 16109b94acceSBarry Smith { 161184cb2905SBarry Smith int ierr; 16129b94acceSBarry Smith int (*r)(SNES); 161377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16147d1a2b2bSBarry Smith if (snes->type == (int) method) return 0; 16157d1a2b2bSBarry Smith 16169b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 161784cb2905SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 1618e3372554SBarry Smith if (!__SNESList) {SETERRQ(1,0,"Could not get methods");} 161937fcc0dbSBarry Smith r = (int (*)(SNES))NRFindRoutine( __SNESList, (int)method, (char *)0 ); 1620e3372554SBarry Smith if (!r) {SETERRQ(1,0,"Unknown method");} 16210452661fSBarry Smith if (snes->data) PetscFree(snes->data); 16229b94acceSBarry Smith snes->set_method_called = 1; 16239b94acceSBarry Smith return (*r)(snes); 16249b94acceSBarry Smith } 16259b94acceSBarry Smith 16269b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16275615d1e5SSatish Balay #undef __FUNC__ 16285eea60f9SBarry Smith #define __FUNC__ "SNESRegister" /* ADIC Ignore */ 16299b94acceSBarry Smith /*@C 16309b94acceSBarry Smith SNESRegister - Adds the method to the nonlinear solver package, given 16314b0e389bSBarry Smith a function pointer and a nonlinear solver name of the type SNESType. 16329b94acceSBarry Smith 16339b94acceSBarry Smith Input Parameters: 16342d872ea7SLois Curfman McInnes . name - either a predefined name such as SNES_EQ_LS, or SNES_NEW 16352d872ea7SLois Curfman McInnes to indicate a new user-defined solver 163640191667SLois Curfman McInnes . sname - corresponding string for name 16379b94acceSBarry Smith . create - routine to create method context 16389b94acceSBarry Smith 163984cb2905SBarry Smith Output Parameter: 164084cb2905SBarry Smith . oname - type associated with this new solver 164184cb2905SBarry Smith 16422d872ea7SLois Curfman McInnes Notes: 16432d872ea7SLois Curfman McInnes Multiple user-defined nonlinear solvers can be added by calling 16442d872ea7SLois Curfman McInnes SNESRegister() with the input parameter "name" set to be SNES_NEW; 16452d872ea7SLois Curfman McInnes each call will return a unique solver type in the output 16462d872ea7SLois Curfman McInnes parameter "oname". 16472d872ea7SLois Curfman McInnes 16489b94acceSBarry Smith .keywords: SNES, nonlinear, register 16499b94acceSBarry Smith 16509b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterDestroy() 16519b94acceSBarry Smith @*/ 165284cb2905SBarry Smith int SNESRegister(SNESType name,SNESType *oname, char *sname, int (*create)(SNES)) 16539b94acceSBarry Smith { 16549b94acceSBarry Smith int ierr; 165584cb2905SBarry Smith static int numberregistered = 0; 165684cb2905SBarry Smith 1657d252947aSBarry Smith if (name == SNES_NEW) name = (SNESType) ((int) SNES_NEW + numberregistered++); 165884cb2905SBarry Smith 165984cb2905SBarry Smith if (oname) *oname = name; 166037fcc0dbSBarry Smith if (!__SNESList) {ierr = NRCreate(&__SNESList); CHKERRQ(ierr);} 166184cb2905SBarry Smith NRRegister( __SNESList, (int) name, sname, (int (*)(void*))create ); 16629b94acceSBarry Smith return 0; 16639b94acceSBarry Smith } 1664a847f771SSatish Balay 16659b94acceSBarry Smith /* --------------------------------------------------------------------- */ 16665615d1e5SSatish Balay #undef __FUNC__ 16675eea60f9SBarry Smith #define __FUNC__ "SNESRegisterDestroy" /* ADIC Ignore */ 16689b94acceSBarry Smith /*@C 16699b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 16709b94acceSBarry Smith registered by SNESRegister(). 16719b94acceSBarry Smith 16729b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 16739b94acceSBarry Smith 16749b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 16759b94acceSBarry Smith @*/ 16769b94acceSBarry Smith int SNESRegisterDestroy() 16779b94acceSBarry Smith { 167837fcc0dbSBarry Smith if (__SNESList) { 167937fcc0dbSBarry Smith NRDestroy( __SNESList ); 168037fcc0dbSBarry Smith __SNESList = 0; 16819b94acceSBarry Smith } 168284cb2905SBarry Smith SNESRegisterAllCalled = 0; 16839b94acceSBarry Smith return 0; 16849b94acceSBarry Smith } 16859b94acceSBarry Smith 16865615d1e5SSatish Balay #undef __FUNC__ 16875eea60f9SBarry Smith #define __FUNC__ "SNESGetTypeFromOptions_Private" /* ADIC Ignore */ 16889b94acceSBarry Smith /* 16894b0e389bSBarry Smith SNESGetTypeFromOptions_Private - Sets the selected method from the 16909b94acceSBarry Smith options database. 16919b94acceSBarry Smith 16929b94acceSBarry Smith Input Parameter: 16939b94acceSBarry Smith . ctx - the SNES context 16949b94acceSBarry Smith 16959b94acceSBarry Smith Output Parameter: 16969b94acceSBarry Smith . method - solver method 16979b94acceSBarry Smith 16989b94acceSBarry Smith Returns: 16999b94acceSBarry Smith Returns 1 if the method is found; 0 otherwise. 17009b94acceSBarry Smith 17019b94acceSBarry Smith Options Database Key: 17024b0e389bSBarry Smith $ -snes_type method 17039b94acceSBarry Smith */ 1704052efed2SBarry Smith int SNESGetTypeFromOptions_Private(SNES ctx,SNESType *method,int *flg) 17059b94acceSBarry Smith { 1706052efed2SBarry Smith int ierr; 17079b94acceSBarry Smith char sbuf[50]; 17085baf8537SBarry Smith 1709052efed2SBarry Smith ierr = OptionsGetString(ctx->prefix,"-snes_type",sbuf,50,flg);CHKERRQ(ierr); 1710052efed2SBarry Smith if (*flg) { 1711052efed2SBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 171237fcc0dbSBarry Smith *method = (SNESType)NRFindID( __SNESList, sbuf ); 17139b94acceSBarry Smith } 17149b94acceSBarry Smith return 0; 17159b94acceSBarry Smith } 17169b94acceSBarry Smith 17175615d1e5SSatish Balay #undef __FUNC__ 17185eea60f9SBarry Smith #define __FUNC__ "SNESGetType" /* ADIC Ignore */ 17199b94acceSBarry Smith /*@C 17209a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 17219b94acceSBarry Smith 17229b94acceSBarry Smith Input Parameter: 17234b0e389bSBarry Smith . snes - nonlinear solver context 17249b94acceSBarry Smith 17259b94acceSBarry Smith Output Parameter: 17269a28b0a6SLois Curfman McInnes . method - SNES method (or use PETSC_NULL) 17279a28b0a6SLois Curfman McInnes . name - name of SNES method (or use PETSC_NULL) 17289b94acceSBarry Smith 17299b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 17309b94acceSBarry Smith @*/ 17314b0e389bSBarry Smith int SNESGetType(SNES snes, SNESType *method,char **name) 17329b94acceSBarry Smith { 173306a9b0adSLois Curfman McInnes int ierr; 173437fcc0dbSBarry Smith if (!__SNESList) {ierr = SNESRegisterAll(); CHKERRQ(ierr);} 17354b0e389bSBarry Smith if (method) *method = (SNESType) snes->type; 173637fcc0dbSBarry Smith if (name) *name = NRFindName( __SNESList, (int) snes->type ); 17379b94acceSBarry Smith return 0; 17389b94acceSBarry Smith } 17399b94acceSBarry Smith 17409b94acceSBarry Smith #include <stdio.h> 17415615d1e5SSatish Balay #undef __FUNC__ 17425eea60f9SBarry Smith #define __FUNC__ "SNESPrintTypes_Private" /* ADIC Ignore */ 17439b94acceSBarry Smith /* 17444b0e389bSBarry Smith SNESPrintTypes_Private - Prints the SNES methods available from the 17459b94acceSBarry Smith options database. 17469b94acceSBarry Smith 17479b94acceSBarry Smith Input Parameters: 174833455573SSatish Balay . comm - communicator (usually MPI_COMM_WORLD) 17499b94acceSBarry Smith . prefix - prefix (usually "-") 17504b0e389bSBarry Smith . name - the options database name (by default "snes_type") 17519b94acceSBarry Smith */ 175233455573SSatish Balay int SNESPrintTypes_Private(MPI_Comm comm,char* prefix,char *name) 17539b94acceSBarry Smith { 17549b94acceSBarry Smith FuncList *entry; 175537fcc0dbSBarry Smith if (!__SNESList) {SNESRegisterAll();} 175637fcc0dbSBarry Smith entry = __SNESList->head; 175777c4ece6SBarry Smith PetscPrintf(comm," %s%s (one of)",prefix,name); 17589b94acceSBarry Smith while (entry) { 175977c4ece6SBarry Smith PetscPrintf(comm," %s",entry->name); 17609b94acceSBarry Smith entry = entry->next; 17619b94acceSBarry Smith } 176277c4ece6SBarry Smith PetscPrintf(comm,"\n"); 17639b94acceSBarry Smith return 0; 17649b94acceSBarry Smith } 17659b94acceSBarry Smith 17665615d1e5SSatish Balay #undef __FUNC__ 17675eea60f9SBarry Smith #define __FUNC__ "SNESGetSolution" /* ADIC Ignore */ 17689b94acceSBarry Smith /*@C 17699b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 17709b94acceSBarry Smith stored. 17719b94acceSBarry Smith 17729b94acceSBarry Smith Input Parameter: 17739b94acceSBarry Smith . snes - the SNES context 17749b94acceSBarry Smith 17759b94acceSBarry Smith Output Parameter: 17769b94acceSBarry Smith . x - the solution 17779b94acceSBarry Smith 17789b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 17799b94acceSBarry Smith 1780a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 17819b94acceSBarry Smith @*/ 17829b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 17839b94acceSBarry Smith { 178477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17859b94acceSBarry Smith *x = snes->vec_sol_always; 17869b94acceSBarry Smith return 0; 17879b94acceSBarry Smith } 17889b94acceSBarry Smith 17895615d1e5SSatish Balay #undef __FUNC__ 17905eea60f9SBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" /* ADIC Ignore */ 17919b94acceSBarry Smith /*@C 17929b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 17939b94acceSBarry Smith stored. 17949b94acceSBarry Smith 17959b94acceSBarry Smith Input Parameter: 17969b94acceSBarry Smith . snes - the SNES context 17979b94acceSBarry Smith 17989b94acceSBarry Smith Output Parameter: 17999b94acceSBarry Smith . x - the solution update 18009b94acceSBarry Smith 18019b94acceSBarry Smith Notes: 18029b94acceSBarry Smith This vector is implementation dependent. 18039b94acceSBarry Smith 18049b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 18059b94acceSBarry Smith 18069b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 18079b94acceSBarry Smith @*/ 18089b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 18099b94acceSBarry Smith { 181077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18119b94acceSBarry Smith *x = snes->vec_sol_update_always; 18129b94acceSBarry Smith return 0; 18139b94acceSBarry Smith } 18149b94acceSBarry Smith 18155615d1e5SSatish Balay #undef __FUNC__ 18165eea60f9SBarry Smith #define __FUNC__ "SNESGetFunction" /* ADIC Ignore */ 18179b94acceSBarry Smith /*@C 18183638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 18199b94acceSBarry Smith 18209b94acceSBarry Smith Input Parameter: 18219b94acceSBarry Smith . snes - the SNES context 18229b94acceSBarry Smith 18239b94acceSBarry Smith Output Parameter: 18243638b69dSLois Curfman McInnes . r - the function 18259b94acceSBarry Smith 18269b94acceSBarry Smith Notes: 18279b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 18289b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 18299b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 18309b94acceSBarry Smith 1831a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 18329b94acceSBarry Smith 183331615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 183431615d04SLois Curfman McInnes SNESGetGradient() 18359b94acceSBarry Smith @*/ 18369b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 18379b94acceSBarry Smith { 183877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1839e3372554SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) SETERRQ(1,0, 18407c792091SSatish Balay "For SNES_NONLINEAR_EQUATIONS only"); 18419b94acceSBarry Smith *r = snes->vec_func_always; 18429b94acceSBarry Smith return 0; 18439b94acceSBarry Smith } 18449b94acceSBarry Smith 18455615d1e5SSatish Balay #undef __FUNC__ 18465eea60f9SBarry Smith #define __FUNC__ "SNESGetGradient" /* ADIC Ignore */ 18479b94acceSBarry Smith /*@C 18483638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 18499b94acceSBarry Smith 18509b94acceSBarry Smith Input Parameter: 18519b94acceSBarry Smith . snes - the SNES context 18529b94acceSBarry Smith 18539b94acceSBarry Smith Output Parameter: 18549b94acceSBarry Smith . r - the gradient 18559b94acceSBarry Smith 18569b94acceSBarry Smith Notes: 18579b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 18589b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18599b94acceSBarry Smith SNESGetFunction(). 18609b94acceSBarry Smith 18619b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 18629b94acceSBarry Smith 186331615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 18649b94acceSBarry Smith @*/ 18659b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 18669b94acceSBarry Smith { 186777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1868e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 186963c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 18709b94acceSBarry Smith *r = snes->vec_func_always; 18719b94acceSBarry Smith return 0; 18729b94acceSBarry Smith } 18739b94acceSBarry Smith 18745615d1e5SSatish Balay #undef __FUNC__ 18755eea60f9SBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" /* ADIC Ignore */ 18769b94acceSBarry Smith /*@ 18779b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 18789b94acceSBarry Smith unconstrained minimization problems. 18799b94acceSBarry Smith 18809b94acceSBarry Smith Input Parameter: 18819b94acceSBarry Smith . snes - the SNES context 18829b94acceSBarry Smith 18839b94acceSBarry Smith Output Parameter: 18849b94acceSBarry Smith . r - the function 18859b94acceSBarry Smith 18869b94acceSBarry Smith Notes: 18879b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 18889b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 18899b94acceSBarry Smith SNESGetFunction(). 18909b94acceSBarry Smith 18919b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 18929b94acceSBarry Smith 189331615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 18949b94acceSBarry Smith @*/ 18959b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 18969b94acceSBarry Smith { 189777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 189874679c65SBarry Smith PetscValidScalarPointer(r); 1899e3372554SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) SETERRQ(1,0, 190063c41f6aSSatish Balay "For SNES_UNCONSTRAINED_MINIMIZATION only"); 19019b94acceSBarry Smith *r = snes->fc; 19029b94acceSBarry Smith return 0; 19039b94acceSBarry Smith } 19049b94acceSBarry Smith 19055615d1e5SSatish Balay #undef __FUNC__ 19065eea60f9SBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" /* ADIC Ignore */ 19073c7409f5SSatish Balay /*@C 19083c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1909d850072dSLois Curfman McInnes SNES options in the database. 19103c7409f5SSatish Balay 19113c7409f5SSatish Balay Input Parameter: 19123c7409f5SSatish Balay . snes - the SNES context 19133c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 19143c7409f5SSatish Balay 1915d850072dSLois Curfman McInnes Notes: 1916a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1917a83b1b31SSatish Balay The first character of all runtime options is AUTOMATICALLY the 1918a83b1b31SSatish Balay hyphen. 1919d850072dSLois Curfman McInnes 19203c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1921a86d99e1SLois Curfman McInnes 1922a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19233c7409f5SSatish Balay @*/ 19243c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 19253c7409f5SSatish Balay { 19263c7409f5SSatish Balay int ierr; 19273c7409f5SSatish Balay 192877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1929639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19303c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19313c7409f5SSatish Balay return 0; 19323c7409f5SSatish Balay } 19333c7409f5SSatish Balay 19345615d1e5SSatish Balay #undef __FUNC__ 19355eea60f9SBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" /* ADIC Ignore */ 19363c7409f5SSatish Balay /*@C 1937f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1938d850072dSLois Curfman McInnes SNES options in the database. 19393c7409f5SSatish Balay 19403c7409f5SSatish Balay Input Parameter: 19413c7409f5SSatish Balay . snes - the SNES context 19423c7409f5SSatish Balay . prefix - the prefix to prepend to all option names 19433c7409f5SSatish Balay 1944d850072dSLois Curfman McInnes Notes: 1945a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1946a83b1b31SSatish Balay The first character of all runtime options is AUTOMATICALLY the 1947a83b1b31SSatish Balay hyphen. 1948d850072dSLois Curfman McInnes 19493c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 1950a86d99e1SLois Curfman McInnes 1951a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 19523c7409f5SSatish Balay @*/ 19533c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 19543c7409f5SSatish Balay { 19553c7409f5SSatish Balay int ierr; 19563c7409f5SSatish Balay 195777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1958639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19593c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19603c7409f5SSatish Balay return 0; 19613c7409f5SSatish Balay } 19623c7409f5SSatish Balay 19635615d1e5SSatish Balay #undef __FUNC__ 19645eea60f9SBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" /* ADIC Ignore */ 1965c83590e2SSatish Balay /*@ 19663c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 19673c7409f5SSatish Balay SNES options in the database. 19683c7409f5SSatish Balay 19693c7409f5SSatish Balay Input Parameter: 19703c7409f5SSatish Balay . snes - the SNES context 19713c7409f5SSatish Balay 19723c7409f5SSatish Balay Output Parameter: 19733c7409f5SSatish Balay . prefix - pointer to the prefix string used 19743c7409f5SSatish Balay 19753c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 1976a86d99e1SLois Curfman McInnes 1977a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 19783c7409f5SSatish Balay @*/ 19793c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 19803c7409f5SSatish Balay { 19813c7409f5SSatish Balay int ierr; 19823c7409f5SSatish Balay 198377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1984639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19853c7409f5SSatish Balay return 0; 19863c7409f5SSatish Balay } 19873c7409f5SSatish Balay 19883c7409f5SSatish Balay 19899b94acceSBarry Smith 19909b94acceSBarry Smith 19919b94acceSBarry Smith 1992