1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*82738288SBarry Smith static char vcid[] = "$Id: snes.c,v 1.159 1998/10/09 19:25:42 bsmith Exp bsmith $"; 39b94acceSBarry Smith #endif 49b94acceSBarry Smith 570f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 6f5eb4b81SSatish Balay #include "src/sys/nreg.h" 79b94acceSBarry Smith #include "pinclude/pviewer.h" 89b94acceSBarry Smith #include <math.h> 99b94acceSBarry Smith 1084cb2905SBarry Smith int SNESRegisterAllCalled = 0; 11488ecbafSBarry Smith FList SNESList = 0; 1282bf6240SBarry Smith 135615d1e5SSatish Balay #undef __FUNC__ 14d4bb536fSBarry Smith #define __FUNC__ "SNESView" 159b94acceSBarry Smith /*@ 169b94acceSBarry Smith SNESView - Prints the SNES data structure. 179b94acceSBarry Smith 18fee21e36SBarry Smith Collective on SNES, unless Viewer is VIEWER_STDOUT_SELF 19fee21e36SBarry Smith 20c7afd0dbSLois Curfman McInnes Input Parameters: 21c7afd0dbSLois Curfman McInnes + SNES - the SNES context 22c7afd0dbSLois Curfman McInnes - viewer - visualization context 23c7afd0dbSLois Curfman McInnes 249b94acceSBarry Smith Options Database Key: 25c8a8ba5cSLois Curfman McInnes . -snes_view - Calls SNESView() at end of SNESSolve() 269b94acceSBarry Smith 279b94acceSBarry Smith Notes: 289b94acceSBarry Smith The available visualization contexts include 29c7afd0dbSLois Curfman McInnes + VIEWER_STDOUT_SELF - standard output (default) 30c7afd0dbSLois Curfman McInnes - VIEWER_STDOUT_WORLD - synchronized standard 31c8a8ba5cSLois Curfman McInnes output where only the first processor opens 32c8a8ba5cSLois Curfman McInnes the file. All other processors send their 33c8a8ba5cSLois Curfman McInnes data to the first processor to print. 349b94acceSBarry Smith 353e081fefSLois Curfman McInnes The user can open an alternative visualization context with 361edb8f11SLois Curfman McInnes 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 513a40ed3dSBarry Smith PetscFunctionBegin; 5274679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5374679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5474679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5574679c65SBarry Smith 5619bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 5719bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 5890ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 5977c4ece6SBarry Smith PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 6082bf6240SBarry Smith SNESGetType(snes,&method); 6177c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," method: %s\n",method); 62e1311b90SBarry Smith if (snes->view) (*snes->view)(snes,viewer); 6377c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 649b94acceSBarry Smith " maximum iterations=%d, maximum function evaluations=%d\n", 659b94acceSBarry Smith snes->max_its,snes->max_funcs); 6677c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 679b94acceSBarry Smith " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 689b94acceSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol); 697a00f4a9SLois Curfman McInnes PetscFPrintf(snes->comm,fd, 707a00f4a9SLois Curfman McInnes " total number of linear solver iterations=%d\n",snes->linear_its); 71ae3c334cSLois Curfman McInnes PetscFPrintf(snes->comm,fd, 7268632166SLois Curfman McInnes " total number of function evaluations=%d\n",snes->nfuncs); 739b94acceSBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 7477c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 759b94acceSBarry Smith if (snes->ksp_ewconv) { 769b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 779b94acceSBarry Smith if (kctx) { 7877c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 799b94acceSBarry Smith " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 809b94acceSBarry Smith kctx->version); 8177c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 829b94acceSBarry Smith " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 839b94acceSBarry Smith kctx->rtol_max,kctx->threshold); 8477c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 859b94acceSBarry Smith kctx->gamma,kctx->alpha,kctx->alpha2); 869b94acceSBarry Smith } 879b94acceSBarry Smith } 88c30f285eSLois Curfman McInnes } else if (vtype == STRING_VIEWER) { 8982bf6240SBarry Smith SNESGetType(snes,&method); 90c30f285eSLois Curfman McInnes ViewerStringSPrintf(viewer," %-3.3s",method); 9119bcc07fSBarry Smith } 929b94acceSBarry Smith SNESGetSLES(snes,&sles); 939b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 943a40ed3dSBarry Smith PetscFunctionReturn(0); 959b94acceSBarry Smith } 969b94acceSBarry Smith 97639f9d9dSBarry Smith /* 98639f9d9dSBarry Smith We retain a list of functions that also take SNES command 99639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 100639f9d9dSBarry Smith */ 101639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 102639f9d9dSBarry Smith static int numberofsetfromoptions; 103639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 104639f9d9dSBarry Smith 1055615d1e5SSatish Balay #undef __FUNC__ 106d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker" 107639f9d9dSBarry Smith /*@ 108639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 109639f9d9dSBarry Smith 110c7afd0dbSLois Curfman McInnes Not Collective 111c7afd0dbSLois Curfman McInnes 112639f9d9dSBarry Smith Input Parameter: 113639f9d9dSBarry Smith . snescheck - function that checks for options 114639f9d9dSBarry Smith 115639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 116639f9d9dSBarry Smith @*/ 117639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 118639f9d9dSBarry Smith { 1193a40ed3dSBarry Smith PetscFunctionBegin; 120639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 121a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed"); 122639f9d9dSBarry Smith } 123639f9d9dSBarry Smith 124639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 1253a40ed3dSBarry Smith PetscFunctionReturn(0); 126639f9d9dSBarry Smith } 127639f9d9dSBarry Smith 1285615d1e5SSatish Balay #undef __FUNC__ 1295615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1309b94acceSBarry Smith /*@ 131682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1329b94acceSBarry Smith 133c7afd0dbSLois Curfman McInnes Collective on SNES 134c7afd0dbSLois Curfman McInnes 1359b94acceSBarry Smith Input Parameter: 1369b94acceSBarry Smith . snes - the SNES context 1379b94acceSBarry Smith 138*82738288SBarry Smith Options Database: 139*82738288SBarry Smith + -snes_type type - SNES_EQ_LS, SNES_EQ_TR, SNES_UM_TR, SNES_UM_LS etc 140*82738288SBarry Smith . -snes_stol - convergence tolerance in terms of the norm 141*82738288SBarry Smith of the change in the solution between steps 142*82738288SBarry Smith . -snes_atol atol - absolute tolerance of residual norm 143*82738288SBarry Smith . -snes_rtol rtol - relative decrease in tolerance norm from initial 144*82738288SBarry Smith . -snes_max_it max_it - maximum number of iterations 145*82738288SBarry Smith . -snes_max_funcs max_funcs - maximum number of function evaluations 146*82738288SBarry Smith . -snes_trtol trtol - trust region tolerance 147*82738288SBarry Smith . -snes_no_convergence_test - skip convergence test in nonlinear or minimization 148*82738288SBarry Smith solver; hence iterations will continue until max_it 149*82738288SBarry Smith or some other criteria is reached. Saves expense 150*82738288SBarry Smith of convergence test 151*82738288SBarry Smith . -snes_monitor - prints residual norm at each iteration 152*82738288SBarry Smith . -snes_xmonitor - plots residual norm at each iteration 153*82738288SBarry Smith - -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 154*82738288SBarry Smith 155*82738288SBarry Smith Options Database for Eisenstat-Walker method: 156*82738288SBarry Smith + -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence 157*82738288SBarry Smith . -snes_ksp_eq_version ver - version of Eisenstat-Walker method 158*82738288SBarry Smith . -snes_ksp_ew_rtol0 rtol0 - 159*82738288SBarry Smith . -snes_ksp_ew_rtolmax rtolmax - 160*82738288SBarry Smith . -snes_ksp_ew_gamma gamma - 161*82738288SBarry Smith . -snes_ksp_ew_alpha alpha - 162*82738288SBarry Smith . -snes_ksp_ew_alpha2 alpha2 - 163*82738288SBarry Smith - -snes_ksp_ew_threshold threshold - 164*82738288SBarry Smith 16511ca99fdSLois Curfman McInnes Notes: 16611ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 16711ca99fdSLois Curfman McInnes the users manual. 16883e2fdc7SBarry Smith 1699b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1709b94acceSBarry Smith 171a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1729b94acceSBarry Smith @*/ 1739b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1749b94acceSBarry Smith { 17582bf6240SBarry Smith char method[256]; 1769b94acceSBarry Smith double tmp; 1779b94acceSBarry Smith SLES sles; 17881f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1793c7409f5SSatish Balay int version = PETSC_DEFAULT; 1809b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1819b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1829b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1839b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1849b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1859b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1869b94acceSBarry Smith 1873a40ed3dSBarry Smith PetscFunctionBegin; 18881f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 18981f57ec7SBarry Smith 19077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 19182bf6240SBarry Smith if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp"); 192ca161407SBarry Smith 19382bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 19482bf6240SBarry Smith ierr = OptionsGetString(snes->prefix,"-snes_type",method,256,&flg); 195052efed2SBarry Smith if (flg) { 19682bf6240SBarry Smith ierr = SNESSetType(snes,(SNESType) method); CHKERRQ(ierr); 1975334005bSBarry Smith } 1983c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 199d64ed03dSBarry Smith if (flg) { 200d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 201d64ed03dSBarry Smith } 2023c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 203d64ed03dSBarry Smith if (flg) { 204d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 205d64ed03dSBarry Smith } 2063c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 207d64ed03dSBarry Smith if (flg) { 208d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 209d64ed03dSBarry Smith } 2103c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 2113c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 212d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 213d64ed03dSBarry Smith if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); } 214d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 215d64ed03dSBarry Smith if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp); CHKERRQ(ierr);} 2163c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 2173c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 218b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 219b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 220b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 221b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 222b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 223b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 224b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 22593c39befSBarry Smith 22693c39befSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr); 22793c39befSBarry Smith if (flg) {snes->converged = 0;} 22893c39befSBarry Smith 2299b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 2309b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 2315bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 2325cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 2333c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 234639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 2353c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 236d31a3109SLois Curfman McInnes if (flg) { 237639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 238d31a3109SLois Curfman McInnes } 23981f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2403c7409f5SSatish Balay if (flg) { 24117699dbbSLois Curfman McInnes int rank = 0; 242d7e8b826SBarry Smith DrawLG lg; 24317699dbbSLois Curfman McInnes MPI_Initialized(&rank); 24417699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 24517699dbbSLois Curfman McInnes if (!rank) { 24681f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2479b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 248f8c826e1SBarry Smith snes->xmonitor = lg; 2499b94acceSBarry Smith PLogObjectParent(snes,lg); 2509b94acceSBarry Smith } 2519b94acceSBarry Smith } 2523c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2533c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2549b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2559b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 256981c4779SBarry Smith PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 257981c4779SBarry Smith } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 25831615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 25931615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 260d64ed03dSBarry Smith PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2619b94acceSBarry Smith } 262639f9d9dSBarry Smith 263639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 264639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 265639f9d9dSBarry Smith } 266639f9d9dSBarry Smith 2679b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2689b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 26993993e2dSLois Curfman McInnes 27082bf6240SBarry Smith /* 27193993e2dSLois Curfman McInnes Since the private setfromoptions requires the type to have 27293993e2dSLois Curfman McInnes been set already, we make sure a type is set by this time. 27382bf6240SBarry Smith */ 27482bf6240SBarry Smith if (!snes->type_name) { 27582bf6240SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 27682bf6240SBarry Smith ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 27782bf6240SBarry Smith } else { 27882bf6240SBarry Smith ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 27982bf6240SBarry Smith } 28082bf6240SBarry Smith } 28182bf6240SBarry Smith 28282bf6240SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 28382bf6240SBarry Smith if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);} 28482bf6240SBarry Smith 2853a40ed3dSBarry Smith if (snes->setfromoptions) { 2863a40ed3dSBarry Smith ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 2873a40ed3dSBarry Smith } 2883a40ed3dSBarry Smith PetscFunctionReturn(0); 2899b94acceSBarry Smith } 2909b94acceSBarry Smith 291a847f771SSatish Balay 2925615d1e5SSatish Balay #undef __FUNC__ 293d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext" 2949b94acceSBarry Smith /*@ 2959b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 2969b94acceSBarry Smith the nonlinear solvers. 2979b94acceSBarry Smith 298fee21e36SBarry Smith Collective on SNES 299fee21e36SBarry Smith 300c7afd0dbSLois Curfman McInnes Input Parameters: 301c7afd0dbSLois Curfman McInnes + snes - the SNES context 302c7afd0dbSLois Curfman McInnes - usrP - optional user context 303c7afd0dbSLois Curfman McInnes 3049b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3059b94acceSBarry Smith 3069b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3079b94acceSBarry Smith @*/ 3089b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3099b94acceSBarry Smith { 3103a40ed3dSBarry Smith PetscFunctionBegin; 31177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3129b94acceSBarry Smith snes->user = usrP; 3133a40ed3dSBarry Smith PetscFunctionReturn(0); 3149b94acceSBarry Smith } 31574679c65SBarry Smith 3165615d1e5SSatish Balay #undef __FUNC__ 317d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext" 3189b94acceSBarry Smith /*@C 3199b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3209b94acceSBarry Smith nonlinear solvers. 3219b94acceSBarry Smith 322c7afd0dbSLois Curfman McInnes Not Collective 323c7afd0dbSLois Curfman McInnes 3249b94acceSBarry Smith Input Parameter: 3259b94acceSBarry Smith . snes - SNES context 3269b94acceSBarry Smith 3279b94acceSBarry Smith Output Parameter: 3289b94acceSBarry Smith . usrP - user context 3299b94acceSBarry Smith 3309b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3319b94acceSBarry Smith 3329b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3339b94acceSBarry Smith @*/ 3349b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3359b94acceSBarry Smith { 3363a40ed3dSBarry Smith PetscFunctionBegin; 33777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3389b94acceSBarry Smith *usrP = snes->user; 3393a40ed3dSBarry Smith PetscFunctionReturn(0); 3409b94acceSBarry Smith } 34174679c65SBarry Smith 3425615d1e5SSatish Balay #undef __FUNC__ 343d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber" 3449b94acceSBarry Smith /*@ 3459b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3469b94acceSBarry Smith nonlinear solver. 3479b94acceSBarry Smith 348c7afd0dbSLois Curfman McInnes Not Collective 349c7afd0dbSLois Curfman McInnes 3509b94acceSBarry Smith Input Parameter: 3519b94acceSBarry Smith . snes - SNES context 3529b94acceSBarry Smith 3539b94acceSBarry Smith Output Parameter: 3549b94acceSBarry Smith . iter - iteration number 3559b94acceSBarry Smith 3569b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3579b94acceSBarry Smith @*/ 3589b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3599b94acceSBarry Smith { 3603a40ed3dSBarry Smith PetscFunctionBegin; 36177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 36274679c65SBarry Smith PetscValidIntPointer(iter); 3639b94acceSBarry Smith *iter = snes->iter; 3643a40ed3dSBarry Smith PetscFunctionReturn(0); 3659b94acceSBarry Smith } 36674679c65SBarry Smith 3675615d1e5SSatish Balay #undef __FUNC__ 3685615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 3699b94acceSBarry Smith /*@ 3709b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3719b94acceSBarry Smith with SNESSSetFunction(). 3729b94acceSBarry Smith 373c7afd0dbSLois Curfman McInnes Collective on SNES 374c7afd0dbSLois Curfman McInnes 3759b94acceSBarry Smith Input Parameter: 3769b94acceSBarry Smith . snes - SNES context 3779b94acceSBarry Smith 3789b94acceSBarry Smith Output Parameter: 3799b94acceSBarry Smith . fnorm - 2-norm of function 3809b94acceSBarry Smith 3819b94acceSBarry Smith Note: 3829b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 383a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 384a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 3859b94acceSBarry Smith 3869b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 387a86d99e1SLois Curfman McInnes 388a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 3899b94acceSBarry Smith @*/ 3909b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 3919b94acceSBarry Smith { 3923a40ed3dSBarry Smith PetscFunctionBegin; 39377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 39474679c65SBarry Smith PetscValidScalarPointer(fnorm); 39574679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 396d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 39774679c65SBarry Smith } 3989b94acceSBarry Smith *fnorm = snes->norm; 3993a40ed3dSBarry Smith PetscFunctionReturn(0); 4009b94acceSBarry Smith } 40174679c65SBarry Smith 4025615d1e5SSatish Balay #undef __FUNC__ 4035615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 4049b94acceSBarry Smith /*@ 4059b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 4069b94acceSBarry Smith with SNESSSetGradient(). 4079b94acceSBarry Smith 408c7afd0dbSLois Curfman McInnes Collective on SNES 409c7afd0dbSLois Curfman McInnes 4109b94acceSBarry Smith Input Parameter: 4119b94acceSBarry Smith . snes - SNES context 4129b94acceSBarry Smith 4139b94acceSBarry Smith Output Parameter: 4149b94acceSBarry Smith . fnorm - 2-norm of gradient 4159b94acceSBarry Smith 4169b94acceSBarry Smith Note: 4179b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 418a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 419a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4209b94acceSBarry Smith 4219b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 422a86d99e1SLois Curfman McInnes 423a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4249b94acceSBarry Smith @*/ 4259b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4269b94acceSBarry Smith { 4273a40ed3dSBarry Smith PetscFunctionBegin; 42877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 42974679c65SBarry Smith PetscValidScalarPointer(gnorm); 43074679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 431d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 43274679c65SBarry Smith } 4339b94acceSBarry Smith *gnorm = snes->norm; 4343a40ed3dSBarry Smith PetscFunctionReturn(0); 4359b94acceSBarry Smith } 43674679c65SBarry Smith 4375615d1e5SSatish Balay #undef __FUNC__ 438d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 4399b94acceSBarry Smith /*@ 4409b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4419b94acceSBarry Smith attempted by the nonlinear solver. 4429b94acceSBarry Smith 443c7afd0dbSLois Curfman McInnes Not Collective 444c7afd0dbSLois Curfman McInnes 4459b94acceSBarry Smith Input Parameter: 4469b94acceSBarry Smith . snes - SNES context 4479b94acceSBarry Smith 4489b94acceSBarry Smith Output Parameter: 4499b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4509b94acceSBarry Smith 451c96a6f78SLois Curfman McInnes Notes: 452c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 453c96a6f78SLois Curfman McInnes 4549b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4559b94acceSBarry Smith @*/ 4569b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4579b94acceSBarry Smith { 4583a40ed3dSBarry Smith PetscFunctionBegin; 45977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 46074679c65SBarry Smith PetscValidIntPointer(nfails); 4619b94acceSBarry Smith *nfails = snes->nfailures; 4623a40ed3dSBarry Smith PetscFunctionReturn(0); 4639b94acceSBarry Smith } 464a847f771SSatish Balay 4655615d1e5SSatish Balay #undef __FUNC__ 466d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" 467c96a6f78SLois Curfman McInnes /*@ 468c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 469c96a6f78SLois Curfman McInnes used by the nonlinear solver. 470c96a6f78SLois Curfman McInnes 471c7afd0dbSLois Curfman McInnes Not Collective 472c7afd0dbSLois Curfman McInnes 473c96a6f78SLois Curfman McInnes Input Parameter: 474c96a6f78SLois Curfman McInnes . snes - SNES context 475c96a6f78SLois Curfman McInnes 476c96a6f78SLois Curfman McInnes Output Parameter: 477c96a6f78SLois Curfman McInnes . lits - number of linear iterations 478c96a6f78SLois Curfman McInnes 479c96a6f78SLois Curfman McInnes Notes: 480c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 481c96a6f78SLois Curfman McInnes 482c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 483c96a6f78SLois Curfman McInnes @*/ 48486bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 485c96a6f78SLois Curfman McInnes { 4863a40ed3dSBarry Smith PetscFunctionBegin; 487c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 488c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 489c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 4903a40ed3dSBarry Smith PetscFunctionReturn(0); 491c96a6f78SLois Curfman McInnes } 492c96a6f78SLois Curfman McInnes 493c96a6f78SLois Curfman McInnes #undef __FUNC__ 494d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES" 4959b94acceSBarry Smith /*@C 4969b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 4979b94acceSBarry Smith 498c7afd0dbSLois Curfman McInnes Not Collective, but if SNES object is parallel, then SLES object is parallel 499c7afd0dbSLois Curfman McInnes 5009b94acceSBarry Smith Input Parameter: 5019b94acceSBarry Smith . snes - the SNES context 5029b94acceSBarry Smith 5039b94acceSBarry Smith Output Parameter: 5049b94acceSBarry Smith . sles - the SLES context 5059b94acceSBarry Smith 5069b94acceSBarry Smith Notes: 5079b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 5089b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5099b94acceSBarry Smith KSP and PC contexts as well. 5109b94acceSBarry Smith 5119b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 5129b94acceSBarry Smith 5139b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 5149b94acceSBarry Smith @*/ 5159b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 5169b94acceSBarry Smith { 5173a40ed3dSBarry Smith PetscFunctionBegin; 51877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5199b94acceSBarry Smith *sles = snes->sles; 5203a40ed3dSBarry Smith PetscFunctionReturn(0); 5219b94acceSBarry Smith } 52282bf6240SBarry Smith 5239b94acceSBarry Smith /* -----------------------------------------------------------*/ 5245615d1e5SSatish Balay #undef __FUNC__ 5255615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 5269b94acceSBarry Smith /*@C 5279b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5289b94acceSBarry Smith 529c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 530c7afd0dbSLois Curfman McInnes 531c7afd0dbSLois Curfman McInnes Input Parameters: 532c7afd0dbSLois Curfman McInnes + comm - MPI communicator 533c7afd0dbSLois Curfman McInnes - type - type of method, either 534c7afd0dbSLois Curfman McInnes SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations) 535c7afd0dbSLois Curfman McInnes or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization) 5369b94acceSBarry Smith 5379b94acceSBarry Smith Output Parameter: 5389b94acceSBarry Smith . outsnes - the new SNES context 5399b94acceSBarry Smith 540c7afd0dbSLois Curfman McInnes Options Database Keys: 541c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 542c7afd0dbSLois Curfman McInnes and no preconditioning matrix 543c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 544c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 545c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 546c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 547c1f60f51SBarry Smith 5489b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5499b94acceSBarry Smith 55063a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 5519b94acceSBarry Smith @*/ 5524b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 5539b94acceSBarry Smith { 5549b94acceSBarry Smith int ierr; 5559b94acceSBarry Smith SNES snes; 5569b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 55737fcc0dbSBarry Smith 5583a40ed3dSBarry Smith PetscFunctionBegin; 5599b94acceSBarry Smith *outsnes = 0; 560d64ed03dSBarry Smith if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){ 561d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 562d64ed03dSBarry Smith } 563f830108cSBarry Smith PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,comm,SNESDestroy,SNESView); 5649b94acceSBarry Smith PLogObjectCreate(snes); 5659b94acceSBarry Smith snes->max_its = 50; 5669750a799SBarry Smith snes->max_funcs = 10000; 5679b94acceSBarry Smith snes->norm = 0.0; 5685a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 569b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 570b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5719b94acceSBarry Smith snes->atol = 1.e-10; 572d64ed03dSBarry Smith } else { 573b4874afaSBarry Smith snes->rtol = 1.e-8; 574b4874afaSBarry Smith snes->ttol = 0.0; 575b4874afaSBarry Smith snes->atol = 1.e-50; 576b4874afaSBarry Smith } 5779b94acceSBarry Smith snes->xtol = 1.e-8; 578b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5799b94acceSBarry Smith snes->nfuncs = 0; 5809b94acceSBarry Smith snes->nfailures = 0; 5817a00f4a9SLois Curfman McInnes snes->linear_its = 0; 582639f9d9dSBarry Smith snes->numbermonitors = 0; 5839b94acceSBarry Smith snes->data = 0; 5849b94acceSBarry Smith snes->view = 0; 5859b94acceSBarry Smith snes->computeumfunction = 0; 5869b94acceSBarry Smith snes->umfunP = 0; 5879b94acceSBarry Smith snes->fc = 0; 5889b94acceSBarry Smith snes->deltatol = 1.e-12; 5899b94acceSBarry Smith snes->fmin = -1.e30; 5909b94acceSBarry Smith snes->method_class = type; 5919b94acceSBarry Smith snes->set_method_called = 0; 59282bf6240SBarry Smith snes->setupcalled = 0; 5939b94acceSBarry Smith snes->ksp_ewconv = 0; 5946f24a144SLois Curfman McInnes snes->xmonitor = 0; 5956f24a144SLois Curfman McInnes snes->vwork = 0; 5966f24a144SLois Curfman McInnes snes->nwork = 0; 597c9005455SLois Curfman McInnes snes->conv_hist_size = 0; 598c9005455SLois Curfman McInnes snes->conv_act_size = 0; 599c9005455SLois Curfman McInnes snes->conv_hist = 0; 6009b94acceSBarry Smith 6019b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 6020452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 603eed86810SBarry Smith PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 6049b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6059b94acceSBarry Smith kctx->version = 2; 6069b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 6079b94acceSBarry Smith this was too large for some test cases */ 6089b94acceSBarry Smith kctx->rtol_last = 0; 6099b94acceSBarry Smith kctx->rtol_max = .9; 6109b94acceSBarry Smith kctx->gamma = 1.0; 6119b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 6129b94acceSBarry Smith kctx->alpha = kctx->alpha2; 6139b94acceSBarry Smith kctx->threshold = .1; 6149b94acceSBarry Smith kctx->lresid_last = 0; 6159b94acceSBarry Smith kctx->norm_last = 0; 6169b94acceSBarry Smith 6179b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 6189b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 6195334005bSBarry Smith 6209b94acceSBarry Smith *outsnes = snes; 6213a40ed3dSBarry Smith PetscFunctionReturn(0); 6229b94acceSBarry Smith } 6239b94acceSBarry Smith 6249b94acceSBarry Smith /* --------------------------------------------------------------- */ 6255615d1e5SSatish Balay #undef __FUNC__ 626d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction" 6279b94acceSBarry Smith /*@C 6289b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6299b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6309b94acceSBarry Smith equations. 6319b94acceSBarry Smith 632fee21e36SBarry Smith Collective on SNES 633fee21e36SBarry Smith 634c7afd0dbSLois Curfman McInnes Input Parameters: 635c7afd0dbSLois Curfman McInnes + snes - the SNES context 636c7afd0dbSLois Curfman McInnes . func - function evaluation routine 637c7afd0dbSLois Curfman McInnes . r - vector to store function value 638c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 639c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6409b94acceSBarry Smith 641c7afd0dbSLois Curfman McInnes Calling sequence of func: 6428d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 643c7afd0dbSLois Curfman McInnes 644313e4042SLois Curfman McInnes . f - function vector 645c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 6469b94acceSBarry Smith 6479b94acceSBarry Smith Notes: 6489b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6499b94acceSBarry Smith $ f'(x) x = -f(x), 650c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 6519b94acceSBarry Smith 6529b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6539b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6549b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 6559b94acceSBarry Smith 6569b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6579b94acceSBarry Smith 658a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6599b94acceSBarry Smith @*/ 6605334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 6619b94acceSBarry Smith { 6623a40ed3dSBarry Smith PetscFunctionBegin; 66377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 664a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 665a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 666a8c6a408SBarry Smith } 6679b94acceSBarry Smith snes->computefunction = func; 6689b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6699b94acceSBarry Smith snes->funP = ctx; 6703a40ed3dSBarry Smith PetscFunctionReturn(0); 6719b94acceSBarry Smith } 6729b94acceSBarry Smith 6735615d1e5SSatish Balay #undef __FUNC__ 6745615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 6759b94acceSBarry Smith /*@ 6769b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6779b94acceSBarry Smith SNESSetFunction(). 6789b94acceSBarry Smith 679c7afd0dbSLois Curfman McInnes Collective on SNES 680c7afd0dbSLois Curfman McInnes 6819b94acceSBarry Smith Input Parameters: 682c7afd0dbSLois Curfman McInnes + snes - the SNES context 683c7afd0dbSLois Curfman McInnes - x - input vector 6849b94acceSBarry Smith 6859b94acceSBarry Smith Output Parameter: 6863638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6879b94acceSBarry Smith 6881bffabb2SLois Curfman McInnes Notes: 6899b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6909b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6919b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6929b94acceSBarry Smith 6939b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6949b94acceSBarry Smith 695a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6969b94acceSBarry Smith @*/ 6979b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 6989b94acceSBarry Smith { 6999b94acceSBarry Smith int ierr; 7009b94acceSBarry Smith 7013a40ed3dSBarry Smith PetscFunctionBegin; 702d4bb536fSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 703a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 704d4bb536fSBarry Smith } 7059b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 706d64ed03dSBarry Smith PetscStackPush("SNES user function"); 7079b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 708d64ed03dSBarry Smith PetscStackPop; 709ae3c334cSLois Curfman McInnes snes->nfuncs++; 7109b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 7113a40ed3dSBarry Smith PetscFunctionReturn(0); 7129b94acceSBarry Smith } 7139b94acceSBarry Smith 7145615d1e5SSatish Balay #undef __FUNC__ 715d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" 7169b94acceSBarry Smith /*@C 7179b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 7189b94acceSBarry Smith unconstrained minimization. 7199b94acceSBarry Smith 720fee21e36SBarry Smith Collective on SNES 721fee21e36SBarry Smith 722c7afd0dbSLois Curfman McInnes Input Parameters: 723c7afd0dbSLois Curfman McInnes + snes - the SNES context 724c7afd0dbSLois Curfman McInnes . func - function evaluation routine 725c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 726c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 7279b94acceSBarry Smith 728c7afd0dbSLois Curfman McInnes Calling sequence of func: 7298d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,double *f,void *ctx); 730c7afd0dbSLois Curfman McInnes 731c7afd0dbSLois Curfman McInnes + x - input vector 7329b94acceSBarry Smith . f - function 733c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined function context 7349b94acceSBarry Smith 7359b94acceSBarry Smith Notes: 7369b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7379b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7389b94acceSBarry Smith SNESSetFunction(). 7399b94acceSBarry Smith 7409b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 7419b94acceSBarry Smith 742a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 743a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 7449b94acceSBarry Smith @*/ 7459b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 7469b94acceSBarry Smith void *ctx) 7479b94acceSBarry Smith { 7483a40ed3dSBarry Smith PetscFunctionBegin; 74977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 750a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 751a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 752a8c6a408SBarry Smith } 7539b94acceSBarry Smith snes->computeumfunction = func; 7549b94acceSBarry Smith snes->umfunP = ctx; 7553a40ed3dSBarry Smith PetscFunctionReturn(0); 7569b94acceSBarry Smith } 7579b94acceSBarry Smith 7585615d1e5SSatish Balay #undef __FUNC__ 7595615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 7609b94acceSBarry Smith /*@ 7619b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 7629b94acceSBarry Smith set with SNESSetMinimizationFunction(). 7639b94acceSBarry Smith 764c7afd0dbSLois Curfman McInnes Collective on SNES 765c7afd0dbSLois Curfman McInnes 7669b94acceSBarry Smith Input Parameters: 767c7afd0dbSLois Curfman McInnes + snes - the SNES context 768c7afd0dbSLois Curfman McInnes - x - input vector 7699b94acceSBarry Smith 7709b94acceSBarry Smith Output Parameter: 7719b94acceSBarry Smith . y - function value 7729b94acceSBarry Smith 7739b94acceSBarry Smith Notes: 7749b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 7759b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7769b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 777a86d99e1SLois Curfman McInnes 778a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 779a86d99e1SLois Curfman McInnes 780a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 781a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 7829b94acceSBarry Smith @*/ 7839b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 7849b94acceSBarry Smith { 7859b94acceSBarry Smith int ierr; 7863a40ed3dSBarry Smith 7873a40ed3dSBarry Smith PetscFunctionBegin; 788a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 789a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 790a8c6a408SBarry Smith } 7919b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 792d64ed03dSBarry Smith PetscStackPush("SNES user minimzation function"); 7939b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 794d64ed03dSBarry Smith PetscStackPop; 795ae3c334cSLois Curfman McInnes snes->nfuncs++; 7969b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7973a40ed3dSBarry Smith PetscFunctionReturn(0); 7989b94acceSBarry Smith } 7999b94acceSBarry Smith 8005615d1e5SSatish Balay #undef __FUNC__ 801d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient" 8029b94acceSBarry Smith /*@C 8039b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 8049b94acceSBarry Smith vector for use by the SNES routines. 8059b94acceSBarry Smith 806c7afd0dbSLois Curfman McInnes Collective on SNES 807c7afd0dbSLois Curfman McInnes 8089b94acceSBarry Smith Input Parameters: 809c7afd0dbSLois Curfman McInnes + snes - the SNES context 8109b94acceSBarry Smith . func - function evaluation routine 811044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 812044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 813c7afd0dbSLois Curfman McInnes - r - vector to store gradient value 814fee21e36SBarry Smith 8159b94acceSBarry Smith Calling sequence of func: 8168d76a1e5SLois Curfman McInnes $ func (SNES, Vec x, Vec g, void *ctx); 8179b94acceSBarry Smith 818c7afd0dbSLois Curfman McInnes + x - input vector 8199b94acceSBarry Smith . g - gradient vector 820c7afd0dbSLois Curfman McInnes - ctx - optional user-defined gradient context 8219b94acceSBarry Smith 8229b94acceSBarry Smith Notes: 8239b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 8249b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 8259b94acceSBarry Smith SNESSetFunction(). 8269b94acceSBarry Smith 8279b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 8289b94acceSBarry Smith 829a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 830a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 8319b94acceSBarry Smith @*/ 83274679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 8339b94acceSBarry Smith { 8343a40ed3dSBarry Smith PetscFunctionBegin; 83577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 836a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 837a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 838a8c6a408SBarry Smith } 8399b94acceSBarry Smith snes->computefunction = func; 8409b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 8419b94acceSBarry Smith snes->funP = ctx; 8423a40ed3dSBarry Smith PetscFunctionReturn(0); 8439b94acceSBarry Smith } 8449b94acceSBarry Smith 8455615d1e5SSatish Balay #undef __FUNC__ 8465615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 8479b94acceSBarry Smith /*@ 848a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 849a86d99e1SLois Curfman McInnes SNESSetGradient(). 8509b94acceSBarry Smith 851c7afd0dbSLois Curfman McInnes Collective on SNES 852c7afd0dbSLois Curfman McInnes 8539b94acceSBarry Smith Input Parameters: 854c7afd0dbSLois Curfman McInnes + snes - the SNES context 855c7afd0dbSLois Curfman McInnes - x - input vector 8569b94acceSBarry Smith 8579b94acceSBarry Smith Output Parameter: 8589b94acceSBarry Smith . y - gradient vector 8599b94acceSBarry Smith 8609b94acceSBarry Smith Notes: 8619b94acceSBarry Smith SNESComputeGradient() is valid only for 8629b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 8639b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 864a86d99e1SLois Curfman McInnes 865a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 866a86d99e1SLois Curfman McInnes 867a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 868a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 8699b94acceSBarry Smith @*/ 8709b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 8719b94acceSBarry Smith { 8729b94acceSBarry Smith int ierr; 8733a40ed3dSBarry Smith 8743a40ed3dSBarry Smith PetscFunctionBegin; 8753a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 876a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8773a40ed3dSBarry Smith } 8789b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 879d64ed03dSBarry Smith PetscStackPush("SNES user gradient function"); 8809b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 881d64ed03dSBarry Smith PetscStackPop; 8829b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 8833a40ed3dSBarry Smith PetscFunctionReturn(0); 8849b94acceSBarry Smith } 8859b94acceSBarry Smith 8865615d1e5SSatish Balay #undef __FUNC__ 8875615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 88862fef451SLois Curfman McInnes /*@ 88962fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 89062fef451SLois Curfman McInnes set with SNESSetJacobian(). 89162fef451SLois Curfman McInnes 892c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 893c7afd0dbSLois Curfman McInnes 89462fef451SLois Curfman McInnes Input Parameters: 895c7afd0dbSLois Curfman McInnes + snes - the SNES context 896c7afd0dbSLois Curfman McInnes - x - input vector 89762fef451SLois Curfman McInnes 89862fef451SLois Curfman McInnes Output Parameters: 899c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 90062fef451SLois Curfman McInnes . B - optional preconditioning matrix 901c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 902fee21e36SBarry Smith 90362fef451SLois Curfman McInnes Notes: 90462fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 90562fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 90662fef451SLois Curfman McInnes 907dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 908dc5a77f8SLois Curfman McInnes flag parameter. 90962fef451SLois Curfman McInnes 91062fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 91162fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 912005c665bSBarry Smith methods is SNESComputeHessian(). 91362fef451SLois Curfman McInnes 91462fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 91562fef451SLois Curfman McInnes 91662fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 91762fef451SLois Curfman McInnes @*/ 9189b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 9199b94acceSBarry Smith { 9209b94acceSBarry Smith int ierr; 9213a40ed3dSBarry Smith 9223a40ed3dSBarry Smith PetscFunctionBegin; 9233a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 924a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 9253a40ed3dSBarry Smith } 9263a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 9279b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 928c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 929d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 9309b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 931d64ed03dSBarry Smith PetscStackPop; 9329b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 9336d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 93477c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 93577c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9363a40ed3dSBarry Smith PetscFunctionReturn(0); 9379b94acceSBarry Smith } 9389b94acceSBarry Smith 9395615d1e5SSatish Balay #undef __FUNC__ 9405615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 94162fef451SLois Curfman McInnes /*@ 94262fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 94362fef451SLois Curfman McInnes set with SNESSetHessian(). 94462fef451SLois Curfman McInnes 945c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 946c7afd0dbSLois Curfman McInnes 94762fef451SLois Curfman McInnes Input Parameters: 948c7afd0dbSLois Curfman McInnes + snes - the SNES context 949c7afd0dbSLois Curfman McInnes - x - input vector 95062fef451SLois Curfman McInnes 95162fef451SLois Curfman McInnes Output Parameters: 952c7afd0dbSLois Curfman McInnes + A - Hessian matrix 95362fef451SLois Curfman McInnes . B - optional preconditioning matrix 954c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 955fee21e36SBarry Smith 95662fef451SLois Curfman McInnes Notes: 95762fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 95862fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 95962fef451SLois Curfman McInnes 960dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 961dc5a77f8SLois Curfman McInnes flag parameter. 96262fef451SLois Curfman McInnes 96362fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 96462fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 96562fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 96662fef451SLois Curfman McInnes 96762fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 96862fef451SLois Curfman McInnes 969a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 970a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 97162fef451SLois Curfman McInnes @*/ 97262fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 9739b94acceSBarry Smith { 9749b94acceSBarry Smith int ierr; 9753a40ed3dSBarry Smith 9763a40ed3dSBarry Smith PetscFunctionBegin; 9773a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 978a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 9793a40ed3dSBarry Smith } 9803a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 98162fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 9820f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 983d64ed03dSBarry Smith PetscStackPush("SNES user Hessian function"); 98462fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 985d64ed03dSBarry Smith PetscStackPop; 98662fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 9870f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 98877c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 98977c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9903a40ed3dSBarry Smith PetscFunctionReturn(0); 9919b94acceSBarry Smith } 9929b94acceSBarry Smith 9935615d1e5SSatish Balay #undef __FUNC__ 994d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian" 9959b94acceSBarry Smith /*@C 9969b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 997044dda88SLois Curfman McInnes location to store the matrix. 9989b94acceSBarry Smith 999c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1000c7afd0dbSLois Curfman McInnes 10019b94acceSBarry Smith Input Parameters: 1002c7afd0dbSLois Curfman McInnes + snes - the SNES context 10039b94acceSBarry Smith . A - Jacobian matrix 10049b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 10059b94acceSBarry Smith . func - Jacobian evaluation routine 1006c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 10072cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 10089b94acceSBarry Smith 10099b94acceSBarry Smith Calling sequence of func: 10108d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 10119b94acceSBarry Smith 1012c7afd0dbSLois Curfman McInnes + x - input vector 10139b94acceSBarry Smith . A - Jacobian matrix 10149b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1015ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1016ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 1017c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 10189b94acceSBarry Smith 10199b94acceSBarry Smith Notes: 1020dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10212cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1022ac21db08SLois Curfman McInnes 1023ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 10249b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 10259b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10269b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10279b94acceSBarry Smith throughout the global iterations. 10289b94acceSBarry Smith 10299b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 10309b94acceSBarry Smith 1031ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 10329b94acceSBarry Smith @*/ 10339b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 10349b94acceSBarry Smith MatStructure*,void*),void *ctx) 10359b94acceSBarry Smith { 10363a40ed3dSBarry Smith PetscFunctionBegin; 103777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1038a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1039a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1040a8c6a408SBarry Smith } 10419b94acceSBarry Smith snes->computejacobian = func; 10429b94acceSBarry Smith snes->jacP = ctx; 10439b94acceSBarry Smith snes->jacobian = A; 10449b94acceSBarry Smith snes->jacobian_pre = B; 10453a40ed3dSBarry Smith PetscFunctionReturn(0); 10469b94acceSBarry Smith } 104762fef451SLois Curfman McInnes 10485615d1e5SSatish Balay #undef __FUNC__ 1049d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian" 1050b4fd4287SBarry Smith /*@ 1051b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1052b4fd4287SBarry Smith provided context for evaluating the Jacobian. 1053b4fd4287SBarry Smith 1054c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 1055c7afd0dbSLois Curfman McInnes 1056b4fd4287SBarry Smith Input Parameter: 1057b4fd4287SBarry Smith . snes - the nonlinear solver context 1058b4fd4287SBarry Smith 1059b4fd4287SBarry Smith Output Parameters: 1060c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 1061b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 1062c7afd0dbSLois Curfman McInnes - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1063fee21e36SBarry Smith 1064b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 1065b4fd4287SBarry Smith @*/ 1066b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1067b4fd4287SBarry Smith { 10683a40ed3dSBarry Smith PetscFunctionBegin; 10693a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1070a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 10713a40ed3dSBarry Smith } 1072b4fd4287SBarry Smith if (A) *A = snes->jacobian; 1073b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 1074b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 10753a40ed3dSBarry Smith PetscFunctionReturn(0); 1076b4fd4287SBarry Smith } 1077b4fd4287SBarry Smith 10785615d1e5SSatish Balay #undef __FUNC__ 1079d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian" 10809b94acceSBarry Smith /*@C 10819b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 1082044dda88SLois Curfman McInnes location to store the matrix. 10839b94acceSBarry Smith 1084c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1085c7afd0dbSLois Curfman McInnes 10869b94acceSBarry Smith Input Parameters: 1087c7afd0dbSLois Curfman McInnes + snes - the SNES context 10889b94acceSBarry Smith . A - Hessian matrix 10899b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 10909b94acceSBarry Smith . func - Jacobian evaluation routine 1091c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 1092313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 10939b94acceSBarry Smith 10949b94acceSBarry Smith Calling sequence of func: 10958d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 10969b94acceSBarry Smith 1097c7afd0dbSLois Curfman McInnes + x - input vector 10989b94acceSBarry Smith . A - Hessian matrix 10999b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1100ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1101ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 1102c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Hessian context 11039b94acceSBarry Smith 11049b94acceSBarry Smith Notes: 1105dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 11062cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1107ac21db08SLois Curfman McInnes 11089b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 11099b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 11109b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 11119b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 11129b94acceSBarry Smith throughout the global iterations. 11139b94acceSBarry Smith 11149b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 11159b94acceSBarry Smith 1116ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 11179b94acceSBarry Smith @*/ 11189b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 11199b94acceSBarry Smith MatStructure*,void*),void *ctx) 11209b94acceSBarry Smith { 11213a40ed3dSBarry Smith PetscFunctionBegin; 112277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1123d4bb536fSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1124a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1125d4bb536fSBarry Smith } 11269b94acceSBarry Smith snes->computejacobian = func; 11279b94acceSBarry Smith snes->jacP = ctx; 11289b94acceSBarry Smith snes->jacobian = A; 11299b94acceSBarry Smith snes->jacobian_pre = B; 11303a40ed3dSBarry Smith PetscFunctionReturn(0); 11319b94acceSBarry Smith } 11329b94acceSBarry Smith 11335615d1e5SSatish Balay #undef __FUNC__ 1134d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian" 113562fef451SLois Curfman McInnes /*@ 113662fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 113762fef451SLois Curfman McInnes provided context for evaluating the Hessian. 113862fef451SLois Curfman McInnes 1139c7afd0dbSLois Curfman McInnes Not Collective, but Mat object is parallel if SNES object is parallel 1140c7afd0dbSLois Curfman McInnes 114162fef451SLois Curfman McInnes Input Parameter: 114262fef451SLois Curfman McInnes . snes - the nonlinear solver context 114362fef451SLois Curfman McInnes 114462fef451SLois Curfman McInnes Output Parameters: 1145c7afd0dbSLois Curfman McInnes + A - location to stash Hessian matrix (or PETSC_NULL) 114662fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 1147c7afd0dbSLois Curfman McInnes - ctx - location to stash Hessian ctx (or PETSC_NULL) 1148fee21e36SBarry Smith 114962fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 1150c7afd0dbSLois Curfman McInnes 1151c7afd0dbSLois Curfman McInnes .keywords: SNES, get, Hessian 115262fef451SLois Curfman McInnes @*/ 115362fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 115462fef451SLois Curfman McInnes { 11553a40ed3dSBarry Smith PetscFunctionBegin; 11563a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){ 1157a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 11583a40ed3dSBarry Smith } 115962fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 116062fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 116162fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 11623a40ed3dSBarry Smith PetscFunctionReturn(0); 116362fef451SLois Curfman McInnes } 116462fef451SLois Curfman McInnes 11659b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 11669b94acceSBarry Smith 11675615d1e5SSatish Balay #undef __FUNC__ 11685615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 11699b94acceSBarry Smith /*@ 11709b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1171272ac6f2SLois Curfman McInnes of a nonlinear solver. 11729b94acceSBarry Smith 1173fee21e36SBarry Smith Collective on SNES 1174fee21e36SBarry Smith 1175c7afd0dbSLois Curfman McInnes Input Parameters: 1176c7afd0dbSLois Curfman McInnes + snes - the SNES context 1177c7afd0dbSLois Curfman McInnes - x - the solution vector 1178c7afd0dbSLois Curfman McInnes 1179272ac6f2SLois Curfman McInnes Notes: 1180272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1181272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1182272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1183272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1184272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1185272ac6f2SLois Curfman McInnes 11869b94acceSBarry Smith .keywords: SNES, nonlinear, setup 11879b94acceSBarry Smith 11889b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 11899b94acceSBarry Smith @*/ 11908ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 11919b94acceSBarry Smith { 1192272ac6f2SLois Curfman McInnes int ierr, flg; 11933a40ed3dSBarry Smith 11943a40ed3dSBarry Smith PetscFunctionBegin; 119577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 119677c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 11978ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 11989b94acceSBarry Smith 1199c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1200c1f60f51SBarry Smith /* 1201c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1202dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1203c1f60f51SBarry Smith */ 1204c1f60f51SBarry Smith if (flg) { 1205c1f60f51SBarry Smith Mat J; 1206c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1207c1f60f51SBarry Smith PLogObjectParent(snes,J); 1208c1f60f51SBarry Smith snes->mfshell = J; 1209c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1210c1f60f51SBarry Smith snes->jacobian = J; 1211c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1212d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1213c1f60f51SBarry Smith snes->jacobian = J; 1214c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1215d4bb536fSBarry Smith } else { 1216a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option"); 1217d4bb536fSBarry Smith } 1218c1f60f51SBarry Smith } 1219272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1220c1f60f51SBarry Smith /* 1221dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1222c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1223c1f60f51SBarry Smith */ 122431615d04SLois Curfman McInnes if (flg) { 1225272ac6f2SLois Curfman McInnes Mat J; 1226272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1227272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1228272ac6f2SLois Curfman McInnes snes->mfshell = J; 122983e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 123083e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 123194a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1232d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 123383e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 123494a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1235d4bb536fSBarry Smith } else { 1236a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option"); 1237d4bb536fSBarry Smith } 1238272ac6f2SLois Curfman McInnes } 12399b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1240a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1241a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1242a8c6a408SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first"); 1243a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 1244a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector"); 1245a8c6a408SBarry Smith } 1246a703fe33SLois Curfman McInnes 1247a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 124882bf6240SBarry Smith if (snes->ksp_ewconv && PetscStrcmp(snes->type_name,SNES_EQ_TR)) { 1249a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1250a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1251a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1252a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1253a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1254a703fe33SLois Curfman McInnes } 12559b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1256a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1257a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1258a8c6a408SBarry Smith if (!snes->computeumfunction) { 1259a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first"); 1260a8c6a408SBarry Smith } 1261a8c6a408SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first"); 1262d4bb536fSBarry Smith } else { 1263a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class"); 1264d4bb536fSBarry Smith } 1265a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 126682bf6240SBarry Smith snes->setupcalled = 1; 12673a40ed3dSBarry Smith PetscFunctionReturn(0); 12689b94acceSBarry Smith } 12699b94acceSBarry Smith 12705615d1e5SSatish Balay #undef __FUNC__ 1271d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy" 12729b94acceSBarry Smith /*@C 12739b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 12749b94acceSBarry Smith with SNESCreate(). 12759b94acceSBarry Smith 1276c7afd0dbSLois Curfman McInnes Collective on SNES 1277c7afd0dbSLois Curfman McInnes 12789b94acceSBarry Smith Input Parameter: 12799b94acceSBarry Smith . snes - the SNES context 12809b94acceSBarry Smith 12819b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 12829b94acceSBarry Smith 128363a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 12849b94acceSBarry Smith @*/ 12859b94acceSBarry Smith int SNESDestroy(SNES snes) 12869b94acceSBarry Smith { 12879b94acceSBarry Smith int ierr; 12883a40ed3dSBarry Smith 12893a40ed3dSBarry Smith PetscFunctionBegin; 129077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12913a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1292d4bb536fSBarry Smith 1293e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes); CHKERRQ(ierr);} 12940452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 1295522c5e43SBarry Smith if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);} 12969b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 1297522c5e43SBarry Smith if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);} 1298522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 12999b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 13000452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 13013a40ed3dSBarry Smith PetscFunctionReturn(0); 13029b94acceSBarry Smith } 13039b94acceSBarry Smith 13049b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 13059b94acceSBarry Smith 13065615d1e5SSatish Balay #undef __FUNC__ 13075615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 13089b94acceSBarry Smith /*@ 1309d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 13109b94acceSBarry Smith 1311c7afd0dbSLois Curfman McInnes Collective on SNES 1312c7afd0dbSLois Curfman McInnes 13139b94acceSBarry Smith Input Parameters: 1314c7afd0dbSLois Curfman McInnes + snes - the SNES context 131533174efeSLois Curfman McInnes . atol - absolute convergence tolerance 131633174efeSLois Curfman McInnes . rtol - relative convergence tolerance 131733174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 131833174efeSLois Curfman McInnes of the change in the solution between steps 131933174efeSLois Curfman McInnes . maxit - maximum number of iterations 1320c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1321fee21e36SBarry Smith 132233174efeSLois Curfman McInnes Options Database Keys: 1323c7afd0dbSLois Curfman McInnes + -snes_atol <atol> - Sets atol 1324c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1325c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1326c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1327c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 13289b94acceSBarry Smith 1329d7a720efSLois Curfman McInnes Notes: 13309b94acceSBarry Smith The default maximum number of iterations is 50. 13319b94acceSBarry Smith The default maximum number of function evaluations is 1000. 13329b94acceSBarry Smith 133333174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 13349b94acceSBarry Smith 1335d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 13369b94acceSBarry Smith @*/ 1337d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 13389b94acceSBarry Smith { 13393a40ed3dSBarry Smith PetscFunctionBegin; 134077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1341d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1342d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1343d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1344d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1345d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 13463a40ed3dSBarry Smith PetscFunctionReturn(0); 13479b94acceSBarry Smith } 13489b94acceSBarry Smith 13495615d1e5SSatish Balay #undef __FUNC__ 13505615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 13519b94acceSBarry Smith /*@ 135233174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 135333174efeSLois Curfman McInnes 1354c7afd0dbSLois Curfman McInnes Not Collective 1355c7afd0dbSLois Curfman McInnes 135633174efeSLois Curfman McInnes Input Parameters: 1357c7afd0dbSLois Curfman McInnes + snes - the SNES context 135833174efeSLois Curfman McInnes . atol - absolute convergence tolerance 135933174efeSLois Curfman McInnes . rtol - relative convergence tolerance 136033174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 136133174efeSLois Curfman McInnes of the change in the solution between steps 136233174efeSLois Curfman McInnes . maxit - maximum number of iterations 1363c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1364fee21e36SBarry Smith 136533174efeSLois Curfman McInnes Notes: 136633174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 136733174efeSLois Curfman McInnes 136833174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 136933174efeSLois Curfman McInnes 137033174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 137133174efeSLois Curfman McInnes @*/ 137233174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 137333174efeSLois Curfman McInnes { 13743a40ed3dSBarry Smith PetscFunctionBegin; 137533174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 137633174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 137733174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 137833174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 137933174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 138033174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 13813a40ed3dSBarry Smith PetscFunctionReturn(0); 138233174efeSLois Curfman McInnes } 138333174efeSLois Curfman McInnes 13845615d1e5SSatish Balay #undef __FUNC__ 13855615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 138633174efeSLois Curfman McInnes /*@ 13879b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 13889b94acceSBarry Smith 1389fee21e36SBarry Smith Collective on SNES 1390fee21e36SBarry Smith 1391c7afd0dbSLois Curfman McInnes Input Parameters: 1392c7afd0dbSLois Curfman McInnes + snes - the SNES context 1393c7afd0dbSLois Curfman McInnes - tol - tolerance 1394c7afd0dbSLois Curfman McInnes 13959b94acceSBarry Smith Options Database Key: 1396c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 13979b94acceSBarry Smith 13989b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 13999b94acceSBarry Smith 1400d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 14019b94acceSBarry Smith @*/ 14029b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 14039b94acceSBarry Smith { 14043a40ed3dSBarry Smith PetscFunctionBegin; 140577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 14069b94acceSBarry Smith snes->deltatol = tol; 14073a40ed3dSBarry Smith PetscFunctionReturn(0); 14089b94acceSBarry Smith } 14099b94acceSBarry Smith 14105615d1e5SSatish Balay #undef __FUNC__ 14115615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 14129b94acceSBarry Smith /*@ 141377c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 14149b94acceSBarry Smith for unconstrained minimization solvers. 14159b94acceSBarry Smith 1416fee21e36SBarry Smith Collective on SNES 1417fee21e36SBarry Smith 1418c7afd0dbSLois Curfman McInnes Input Parameters: 1419c7afd0dbSLois Curfman McInnes + snes - the SNES context 1420c7afd0dbSLois Curfman McInnes - ftol - minimum function tolerance 1421c7afd0dbSLois Curfman McInnes 14229b94acceSBarry Smith Options Database Key: 1423c7afd0dbSLois Curfman McInnes . -snes_fmin <ftol> - Sets ftol 14249b94acceSBarry Smith 14259b94acceSBarry Smith Note: 142677c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 14279b94acceSBarry Smith methods only. 14289b94acceSBarry Smith 14299b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 14309b94acceSBarry Smith 1431d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 14329b94acceSBarry Smith @*/ 143377c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 14349b94acceSBarry Smith { 14353a40ed3dSBarry Smith PetscFunctionBegin; 143677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 14379b94acceSBarry Smith snes->fmin = ftol; 14383a40ed3dSBarry Smith PetscFunctionReturn(0); 14399b94acceSBarry Smith } 14409b94acceSBarry Smith 14419b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 14429b94acceSBarry Smith 14435615d1e5SSatish Balay #undef __FUNC__ 1444d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor" 14459b94acceSBarry Smith /*@C 14465cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 14479b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 14489b94acceSBarry Smith progress. 14499b94acceSBarry Smith 1450fee21e36SBarry Smith Collective on SNES 1451fee21e36SBarry Smith 1452c7afd0dbSLois Curfman McInnes Input Parameters: 1453c7afd0dbSLois Curfman McInnes + snes - the SNES context 1454c7afd0dbSLois Curfman McInnes . func - monitoring routine 1455c7afd0dbSLois Curfman McInnes - mctx - [optional] user-defined context for private data for the 1456c7afd0dbSLois Curfman McInnes monitor routine (may be PETSC_NULL) 14579b94acceSBarry Smith 1458c7afd0dbSLois Curfman McInnes Calling sequence of func: 14598d76a1e5SLois Curfman McInnes $ int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 1460c7afd0dbSLois Curfman McInnes 1461c7afd0dbSLois Curfman McInnes + snes - the SNES context 1462c7afd0dbSLois Curfman McInnes . its - iteration number 1463c7afd0dbSLois Curfman McInnes . mctx - [optional] monitoring context 1464c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 1465c7afd0dbSLois Curfman McInnes for SNES_NONLINEAR_EQUATIONS methods 1466c7afd0dbSLois Curfman McInnes - norm - 2-norm gradient value (may be estimated) 1467c7afd0dbSLois Curfman McInnes for SNES_UNCONSTRAINED_MINIMIZATION methods 14689b94acceSBarry Smith 14699665c990SLois Curfman McInnes Options Database Keys: 1470c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1471c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1472c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1473c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1474c7afd0dbSLois Curfman McInnes been hardwired into a code by 1475c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1476c7afd0dbSLois Curfman McInnes does not cancel those set via 1477c7afd0dbSLois Curfman McInnes the options database. 14789665c990SLois Curfman McInnes 1479639f9d9dSBarry Smith Notes: 14806bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 14816bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 14826bc08f3fSLois Curfman McInnes order in which they were set. 1483639f9d9dSBarry Smith 14849b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 14859b94acceSBarry Smith 14865cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 14879b94acceSBarry Smith @*/ 148874679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 14899b94acceSBarry Smith { 14903a40ed3dSBarry Smith PetscFunctionBegin; 1491639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1492a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set"); 1493639f9d9dSBarry Smith } 1494639f9d9dSBarry Smith 1495639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1496639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 14973a40ed3dSBarry Smith PetscFunctionReturn(0); 14989b94acceSBarry Smith } 14999b94acceSBarry Smith 15005615d1e5SSatish Balay #undef __FUNC__ 15015cd90555SBarry Smith #define __FUNC__ "SNESClearMonitor" 15025cd90555SBarry Smith /*@C 15035cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 15045cd90555SBarry Smith 1505c7afd0dbSLois Curfman McInnes Collective on SNES 1506c7afd0dbSLois Curfman McInnes 15075cd90555SBarry Smith Input Parameters: 15085cd90555SBarry Smith . snes - the SNES context 15095cd90555SBarry Smith 15105cd90555SBarry Smith Options Database: 1511c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1512c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1513c7afd0dbSLois Curfman McInnes set via the options database 15145cd90555SBarry Smith 15155cd90555SBarry Smith Notes: 15165cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 15175cd90555SBarry Smith 15185cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 15195cd90555SBarry Smith 15205cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 15215cd90555SBarry Smith @*/ 15225cd90555SBarry Smith int SNESClearMonitor( SNES snes ) 15235cd90555SBarry Smith { 15245cd90555SBarry Smith PetscFunctionBegin; 15255cd90555SBarry Smith snes->numbermonitors = 0; 15265cd90555SBarry Smith PetscFunctionReturn(0); 15275cd90555SBarry Smith } 15285cd90555SBarry Smith 15295cd90555SBarry Smith #undef __FUNC__ 1530d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest" 15319b94acceSBarry Smith /*@C 15329b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 15339b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 15349b94acceSBarry Smith 1535fee21e36SBarry Smith Collective on SNES 1536fee21e36SBarry Smith 1537c7afd0dbSLois Curfman McInnes Input Parameters: 1538c7afd0dbSLois Curfman McInnes + snes - the SNES context 1539c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1540c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1541c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 15429b94acceSBarry Smith 1543c7afd0dbSLois Curfman McInnes Calling sequence of func: 15448d76a1e5SLois Curfman McInnes $ int func (SNES snes,double xnorm,double gnorm,double f,void *cctx) 1545c7afd0dbSLois Curfman McInnes 1546c7afd0dbSLois Curfman McInnes + snes - the SNES context 1547c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1548c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 1549c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods) 1550c7afd0dbSLois Curfman McInnes . f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods) 1551c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods) 1552c7afd0dbSLois Curfman McInnes - f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods) 15539b94acceSBarry Smith 15549b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 15559b94acceSBarry Smith 155640191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 155740191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 15589b94acceSBarry Smith @*/ 155974679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 15609b94acceSBarry Smith { 15613a40ed3dSBarry Smith PetscFunctionBegin; 15629b94acceSBarry Smith (snes)->converged = func; 15639b94acceSBarry Smith (snes)->cnvP = cctx; 15643a40ed3dSBarry Smith PetscFunctionReturn(0); 15659b94acceSBarry Smith } 15669b94acceSBarry Smith 15675615d1e5SSatish Balay #undef __FUNC__ 1568d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" 1569c9005455SLois Curfman McInnes /*@ 1570c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1571c9005455SLois Curfman McInnes 1572fee21e36SBarry Smith Collective on SNES 1573fee21e36SBarry Smith 1574c7afd0dbSLois Curfman McInnes Input Parameters: 1575c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1576c7afd0dbSLois Curfman McInnes . a - array to hold history 1577c7afd0dbSLois Curfman McInnes - na - size of a 1578c7afd0dbSLois Curfman McInnes 1579c9005455SLois Curfman McInnes Notes: 1580c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1581c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1582c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1583c9005455SLois Curfman McInnes at each step. 1584c9005455SLois Curfman McInnes 1585c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1586c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1587c9005455SLois Curfman McInnes during the section of code that is being timed. 1588c9005455SLois Curfman McInnes 1589c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1590c9005455SLois Curfman McInnes @*/ 1591c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na) 1592c9005455SLois Curfman McInnes { 15933a40ed3dSBarry Smith PetscFunctionBegin; 1594c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1595c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1596c9005455SLois Curfman McInnes snes->conv_hist = a; 1597c9005455SLois Curfman McInnes snes->conv_hist_size = na; 15983a40ed3dSBarry Smith PetscFunctionReturn(0); 1599c9005455SLois Curfman McInnes } 1600c9005455SLois Curfman McInnes 1601c9005455SLois Curfman McInnes #undef __FUNC__ 16025615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 16039b94acceSBarry Smith /* 16049b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 16059b94acceSBarry Smith positive parameter delta. 16069b94acceSBarry Smith 16079b94acceSBarry Smith Input Parameters: 1608c7afd0dbSLois Curfman McInnes + snes - the SNES context 16099b94acceSBarry Smith . y - approximate solution of linear system 16109b94acceSBarry Smith . fnorm - 2-norm of current function 1611c7afd0dbSLois Curfman McInnes - delta - trust region size 16129b94acceSBarry Smith 16139b94acceSBarry Smith Output Parameters: 1614c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 16159b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 16169b94acceSBarry Smith region, and exceeds zero otherwise. 1617c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 16189b94acceSBarry Smith 16199b94acceSBarry Smith Note: 162040191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 16219b94acceSBarry Smith is set to be the maximum allowable step size. 16229b94acceSBarry Smith 16239b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16249b94acceSBarry Smith */ 16259b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 16269b94acceSBarry Smith double *gpnorm,double *ynorm) 16279b94acceSBarry Smith { 16289b94acceSBarry Smith double norm; 16299b94acceSBarry Smith Scalar cnorm; 16303a40ed3dSBarry Smith int ierr; 16313a40ed3dSBarry Smith 16323a40ed3dSBarry Smith PetscFunctionBegin; 16333a40ed3dSBarry Smith ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr); 16349b94acceSBarry Smith if (norm > *delta) { 16359b94acceSBarry Smith norm = *delta/norm; 16369b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 16379b94acceSBarry Smith cnorm = norm; 16389b94acceSBarry Smith VecScale( &cnorm, y ); 16399b94acceSBarry Smith *ynorm = *delta; 16409b94acceSBarry Smith } else { 16419b94acceSBarry Smith *gpnorm = 0.0; 16429b94acceSBarry Smith *ynorm = norm; 16439b94acceSBarry Smith } 16443a40ed3dSBarry Smith PetscFunctionReturn(0); 16459b94acceSBarry Smith } 16469b94acceSBarry Smith 16475615d1e5SSatish Balay #undef __FUNC__ 16485615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 16499b94acceSBarry Smith /*@ 16509b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 165163a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 16529b94acceSBarry Smith 1653c7afd0dbSLois Curfman McInnes Collective on SNES 1654c7afd0dbSLois Curfman McInnes 1655b2002411SLois Curfman McInnes Input Parameters: 1656c7afd0dbSLois Curfman McInnes + snes - the SNES context 1657c7afd0dbSLois Curfman McInnes - x - the solution vector 16589b94acceSBarry Smith 16599b94acceSBarry Smith Output Parameter: 1660b2002411SLois Curfman McInnes . its - number of iterations until termination 16619b94acceSBarry Smith 1662b2002411SLois Curfman McInnes Notes: 16638ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 16648ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 16658ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 16668ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 16678ddd3da0SLois Curfman McInnes 16689b94acceSBarry Smith .keywords: SNES, nonlinear, solve 16699b94acceSBarry Smith 167063a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 16719b94acceSBarry Smith @*/ 16728ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 16739b94acceSBarry Smith { 16743c7409f5SSatish Balay int ierr, flg; 1675052efed2SBarry Smith 16763a40ed3dSBarry Smith PetscFunctionBegin; 167777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 167874679c65SBarry Smith PetscValidIntPointer(its); 167982bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1680c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 16819b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 1682c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 16839b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 16849b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 16853c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 16866d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 16873a40ed3dSBarry Smith PetscFunctionReturn(0); 16889b94acceSBarry Smith } 16899b94acceSBarry Smith 16909b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 16919b94acceSBarry Smith 16925615d1e5SSatish Balay #undef __FUNC__ 16935615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 169482bf6240SBarry Smith /*@C 16954b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 16969b94acceSBarry Smith 1697fee21e36SBarry Smith Collective on SNES 1698fee21e36SBarry Smith 1699c7afd0dbSLois Curfman McInnes Input Parameters: 1700c7afd0dbSLois Curfman McInnes + snes - the SNES context 1701c7afd0dbSLois Curfman McInnes - method - a known method 1702c7afd0dbSLois Curfman McInnes 1703c7afd0dbSLois Curfman McInnes Options Database Key: 1704c7afd0dbSLois Curfman McInnes . -snes_type <method> - Sets the method; use -help for a list 1705c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1706ae12b187SLois Curfman McInnes 17079b94acceSBarry Smith Notes: 17089b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 1709c7afd0dbSLois Curfman McInnes . SNES_EQ_LS - Newton's method with line search 1710c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 1711c7afd0dbSLois Curfman McInnes . SNES_EQ_TR - Newton's method with trust region 1712c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 1713c7afd0dbSLois Curfman McInnes . SNES_UM_TR - Newton's method with trust region 1714c7afd0dbSLois Curfman McInnes (unconstrained minimization) 1715c7afd0dbSLois Curfman McInnes . SNES_UM_LS - Newton's method with line search 1716c7afd0dbSLois Curfman McInnes (unconstrained minimization) 17179b94acceSBarry Smith 1718ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1719ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1720ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1721ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1722ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1723ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1724ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1725ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1726ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1727ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1728a703fe33SLois Curfman McInnes 1729f536c699SSatish Balay .keywords: SNES, set, method 17309b94acceSBarry Smith @*/ 17314b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 17329b94acceSBarry Smith { 173384cb2905SBarry Smith int ierr; 17349b94acceSBarry Smith int (*r)(SNES); 17353a40ed3dSBarry Smith 17363a40ed3dSBarry Smith PetscFunctionBegin; 173777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 173882bf6240SBarry Smith 173982bf6240SBarry Smith if (!PetscStrcmp(snes->type_name,method)) PetscFunctionReturn(0); 174082bf6240SBarry Smith 174182bf6240SBarry Smith if (snes->setupcalled) { 1742e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 174382bf6240SBarry Smith snes->data = 0; 174482bf6240SBarry Smith } 17457d1a2b2bSBarry Smith 17469b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 174782bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);} 174882bf6240SBarry Smith 1749488ecbafSBarry Smith ierr = FListFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr); 175082bf6240SBarry Smith 17511548b14aSBarry Smith if (!r) SETERRQ(1,1,"Unable to find requested SNES type"); 17521548b14aSBarry Smith 17530452661fSBarry Smith if (snes->data) PetscFree(snes->data); 175482bf6240SBarry Smith snes->data = 0; 17553a40ed3dSBarry Smith ierr = (*r)(snes); CHKERRQ(ierr); 175682bf6240SBarry Smith 175782bf6240SBarry Smith if (snes->type_name) PetscFree(snes->type_name); 175882bf6240SBarry Smith snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name); 175982bf6240SBarry Smith PetscStrcpy(snes->type_name,method); 176082bf6240SBarry Smith snes->set_method_called = 1; 176182bf6240SBarry Smith 17623a40ed3dSBarry Smith PetscFunctionReturn(0); 17639b94acceSBarry Smith } 17649b94acceSBarry Smith 1765a847f771SSatish Balay 17669b94acceSBarry Smith /* --------------------------------------------------------------------- */ 17675615d1e5SSatish Balay #undef __FUNC__ 1768d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy" 17699b94acceSBarry Smith /*@C 17709b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 17719b94acceSBarry Smith registered by SNESRegister(). 17729b94acceSBarry Smith 1773fee21e36SBarry Smith Not Collective 1774fee21e36SBarry Smith 17759b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 17769b94acceSBarry Smith 17779b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 17789b94acceSBarry Smith @*/ 1779cf256101SBarry Smith int SNESRegisterDestroy(void) 17809b94acceSBarry Smith { 178182bf6240SBarry Smith int ierr; 178282bf6240SBarry Smith 17833a40ed3dSBarry Smith PetscFunctionBegin; 178482bf6240SBarry Smith if (SNESList) { 1785488ecbafSBarry Smith ierr = FListDestroy( SNESList );CHKERRQ(ierr); 178682bf6240SBarry Smith SNESList = 0; 17879b94acceSBarry Smith } 178884cb2905SBarry Smith SNESRegisterAllCalled = 0; 17893a40ed3dSBarry Smith PetscFunctionReturn(0); 17909b94acceSBarry Smith } 17919b94acceSBarry Smith 17925615d1e5SSatish Balay #undef __FUNC__ 1793d4bb536fSBarry Smith #define __FUNC__ "SNESGetType" 17949b94acceSBarry Smith /*@C 17959a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 17969b94acceSBarry Smith 1797c7afd0dbSLois Curfman McInnes Not Collective 1798c7afd0dbSLois Curfman McInnes 17999b94acceSBarry Smith Input Parameter: 18004b0e389bSBarry Smith . snes - nonlinear solver context 18019b94acceSBarry Smith 18029b94acceSBarry Smith Output Parameter: 180382bf6240SBarry Smith . method - SNES method (a charactor string) 18049b94acceSBarry Smith 18059b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 18069b94acceSBarry Smith @*/ 180782bf6240SBarry Smith int SNESGetType(SNES snes, SNESType *method) 18089b94acceSBarry Smith { 18093a40ed3dSBarry Smith PetscFunctionBegin; 181082bf6240SBarry Smith *method = snes->type_name; 18113a40ed3dSBarry Smith PetscFunctionReturn(0); 18129b94acceSBarry Smith } 18139b94acceSBarry Smith 18145615d1e5SSatish Balay #undef __FUNC__ 1815d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution" 18169b94acceSBarry Smith /*@C 18179b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 18189b94acceSBarry Smith stored. 18199b94acceSBarry Smith 1820c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1821c7afd0dbSLois Curfman McInnes 18229b94acceSBarry Smith Input Parameter: 18239b94acceSBarry Smith . snes - the SNES context 18249b94acceSBarry Smith 18259b94acceSBarry Smith Output Parameter: 18269b94acceSBarry Smith . x - the solution 18279b94acceSBarry Smith 18289b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 18299b94acceSBarry Smith 1830a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 18319b94acceSBarry Smith @*/ 18329b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 18339b94acceSBarry Smith { 18343a40ed3dSBarry Smith PetscFunctionBegin; 183577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18369b94acceSBarry Smith *x = snes->vec_sol_always; 18373a40ed3dSBarry Smith PetscFunctionReturn(0); 18389b94acceSBarry Smith } 18399b94acceSBarry Smith 18405615d1e5SSatish Balay #undef __FUNC__ 1841d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" 18429b94acceSBarry Smith /*@C 18439b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 18449b94acceSBarry Smith stored. 18459b94acceSBarry Smith 1846c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1847c7afd0dbSLois Curfman McInnes 18489b94acceSBarry Smith Input Parameter: 18499b94acceSBarry Smith . snes - the SNES context 18509b94acceSBarry Smith 18519b94acceSBarry Smith Output Parameter: 18529b94acceSBarry Smith . x - the solution update 18539b94acceSBarry Smith 18549b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 18559b94acceSBarry Smith 18569b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 18579b94acceSBarry Smith @*/ 18589b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 18599b94acceSBarry Smith { 18603a40ed3dSBarry Smith PetscFunctionBegin; 186177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18629b94acceSBarry Smith *x = snes->vec_sol_update_always; 18633a40ed3dSBarry Smith PetscFunctionReturn(0); 18649b94acceSBarry Smith } 18659b94acceSBarry Smith 18665615d1e5SSatish Balay #undef __FUNC__ 1867d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction" 18689b94acceSBarry Smith /*@C 18693638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 18709b94acceSBarry Smith 1871c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1872c7afd0dbSLois Curfman McInnes 18739b94acceSBarry Smith Input Parameter: 18749b94acceSBarry Smith . snes - the SNES context 18759b94acceSBarry Smith 18769b94acceSBarry Smith Output Parameter: 18773638b69dSLois Curfman McInnes . r - the function 18789b94acceSBarry Smith 18799b94acceSBarry Smith Notes: 18809b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 18819b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 18829b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 18839b94acceSBarry Smith 1884a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 18859b94acceSBarry Smith 188631615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 188731615d04SLois Curfman McInnes SNESGetGradient() 18889b94acceSBarry Smith @*/ 18899b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 18909b94acceSBarry Smith { 18913a40ed3dSBarry Smith PetscFunctionBegin; 189277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1893a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1894a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1895a8c6a408SBarry Smith } 18969b94acceSBarry Smith *r = snes->vec_func_always; 18973a40ed3dSBarry Smith PetscFunctionReturn(0); 18989b94acceSBarry Smith } 18999b94acceSBarry Smith 19005615d1e5SSatish Balay #undef __FUNC__ 1901d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient" 19029b94acceSBarry Smith /*@C 19033638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 19049b94acceSBarry Smith 1905c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1906c7afd0dbSLois Curfman McInnes 19079b94acceSBarry Smith Input Parameter: 19089b94acceSBarry Smith . snes - the SNES context 19099b94acceSBarry Smith 19109b94acceSBarry Smith Output Parameter: 19119b94acceSBarry Smith . r - the gradient 19129b94acceSBarry Smith 19139b94acceSBarry Smith Notes: 19149b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 19159b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 19169b94acceSBarry Smith SNESGetFunction(). 19179b94acceSBarry Smith 19189b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 19199b94acceSBarry Smith 192031615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 19219b94acceSBarry Smith @*/ 19229b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 19239b94acceSBarry Smith { 19243a40ed3dSBarry Smith PetscFunctionBegin; 192577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 19263a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1927a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 19283a40ed3dSBarry Smith } 19299b94acceSBarry Smith *r = snes->vec_func_always; 19303a40ed3dSBarry Smith PetscFunctionReturn(0); 19319b94acceSBarry Smith } 19329b94acceSBarry Smith 19335615d1e5SSatish Balay #undef __FUNC__ 1934d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" 19359b94acceSBarry Smith /*@ 19369b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 19379b94acceSBarry Smith unconstrained minimization problems. 19389b94acceSBarry Smith 1939c7afd0dbSLois Curfman McInnes Not Collective 1940c7afd0dbSLois Curfman McInnes 19419b94acceSBarry Smith Input Parameter: 19429b94acceSBarry Smith . snes - the SNES context 19439b94acceSBarry Smith 19449b94acceSBarry Smith Output Parameter: 19459b94acceSBarry Smith . r - the function 19469b94acceSBarry Smith 19479b94acceSBarry Smith Notes: 19489b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 19499b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 19509b94acceSBarry Smith SNESGetFunction(). 19519b94acceSBarry Smith 19529b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 19539b94acceSBarry Smith 195431615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 19559b94acceSBarry Smith @*/ 19569b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 19579b94acceSBarry Smith { 19583a40ed3dSBarry Smith PetscFunctionBegin; 195977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 196074679c65SBarry Smith PetscValidScalarPointer(r); 19613a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1962a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 19633a40ed3dSBarry Smith } 19649b94acceSBarry Smith *r = snes->fc; 19653a40ed3dSBarry Smith PetscFunctionReturn(0); 19669b94acceSBarry Smith } 19679b94acceSBarry Smith 19685615d1e5SSatish Balay #undef __FUNC__ 1969d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" 19703c7409f5SSatish Balay /*@C 19713c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1972d850072dSLois Curfman McInnes SNES options in the database. 19733c7409f5SSatish Balay 1974fee21e36SBarry Smith Collective on SNES 1975fee21e36SBarry Smith 1976c7afd0dbSLois Curfman McInnes Input Parameter: 1977c7afd0dbSLois Curfman McInnes + snes - the SNES context 1978c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1979c7afd0dbSLois Curfman McInnes 1980d850072dSLois Curfman McInnes Notes: 1981a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1982c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1983d850072dSLois Curfman McInnes 19843c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1985a86d99e1SLois Curfman McInnes 1986a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19873c7409f5SSatish Balay @*/ 19883c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 19893c7409f5SSatish Balay { 19903c7409f5SSatish Balay int ierr; 19913c7409f5SSatish Balay 19923a40ed3dSBarry Smith PetscFunctionBegin; 199377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1994639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19953c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19963a40ed3dSBarry Smith PetscFunctionReturn(0); 19973c7409f5SSatish Balay } 19983c7409f5SSatish Balay 19995615d1e5SSatish Balay #undef __FUNC__ 2000d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" 20013c7409f5SSatish Balay /*@C 2002f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2003d850072dSLois Curfman McInnes SNES options in the database. 20043c7409f5SSatish Balay 2005fee21e36SBarry Smith Collective on SNES 2006fee21e36SBarry Smith 2007c7afd0dbSLois Curfman McInnes Input Parameters: 2008c7afd0dbSLois Curfman McInnes + snes - the SNES context 2009c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2010c7afd0dbSLois Curfman McInnes 2011d850072dSLois Curfman McInnes Notes: 2012a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2013c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2014d850072dSLois Curfman McInnes 20153c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2016a86d99e1SLois Curfman McInnes 2017a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20183c7409f5SSatish Balay @*/ 20193c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 20203c7409f5SSatish Balay { 20213c7409f5SSatish Balay int ierr; 20223c7409f5SSatish Balay 20233a40ed3dSBarry Smith PetscFunctionBegin; 202477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2025639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 20263c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 20273a40ed3dSBarry Smith PetscFunctionReturn(0); 20283c7409f5SSatish Balay } 20293c7409f5SSatish Balay 20305615d1e5SSatish Balay #undef __FUNC__ 2031d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" 2032c83590e2SSatish Balay /*@ 20333c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20343c7409f5SSatish Balay SNES options in the database. 20353c7409f5SSatish Balay 2036c7afd0dbSLois Curfman McInnes Not Collective 2037c7afd0dbSLois Curfman McInnes 20383c7409f5SSatish Balay Input Parameter: 20393c7409f5SSatish Balay . snes - the SNES context 20403c7409f5SSatish Balay 20413c7409f5SSatish Balay Output Parameter: 20423c7409f5SSatish Balay . prefix - pointer to the prefix string used 20433c7409f5SSatish Balay 20443c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2045a86d99e1SLois Curfman McInnes 2046a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20473c7409f5SSatish Balay @*/ 20483c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 20493c7409f5SSatish Balay { 20503c7409f5SSatish Balay int ierr; 20513c7409f5SSatish Balay 20523a40ed3dSBarry Smith PetscFunctionBegin; 205377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2054639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 20553a40ed3dSBarry Smith PetscFunctionReturn(0); 20563c7409f5SSatish Balay } 20573c7409f5SSatish Balay 2058ca161407SBarry Smith #undef __FUNC__ 2059ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp" 2060ca161407SBarry Smith /*@ 2061ca161407SBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2062ca161407SBarry Smith 2063c7afd0dbSLois Curfman McInnes Collective on SNES 2064c7afd0dbSLois Curfman McInnes 2065ca161407SBarry Smith Input Parameter: 2066ca161407SBarry Smith . snes - the SNES context 2067ca161407SBarry Smith 2068ca161407SBarry Smith Options Database Keys: 2069c7afd0dbSLois Curfman McInnes + -help - Prints SNES options 2070c7afd0dbSLois Curfman McInnes - -h - Prints SNES options 2071ca161407SBarry Smith 2072ca161407SBarry Smith .keywords: SNES, nonlinear, help 2073ca161407SBarry Smith 2074ca161407SBarry Smith .seealso: SNESSetFromOptions() 2075ca161407SBarry Smith @*/ 2076ca161407SBarry Smith int SNESPrintHelp(SNES snes) 2077ca161407SBarry Smith { 2078ca161407SBarry Smith char p[64]; 2079ca161407SBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2080ca161407SBarry Smith int ierr; 2081ca161407SBarry Smith 2082ca161407SBarry Smith PetscFunctionBegin; 2083ca161407SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2084ca161407SBarry Smith 2085ca161407SBarry Smith PetscStrcpy(p,"-"); 2086ca161407SBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2087ca161407SBarry Smith 2088ca161407SBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2089ca161407SBarry Smith 209076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n"); 2091488ecbafSBarry Smith ierr = FListPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr); 209276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 209376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 209476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 209576be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 209676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 209776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 209876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 209976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 210076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 210176be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 2102ca161407SBarry Smith residual norm at each iteration.\n",p); 210376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 2104ca161407SBarry Smith residual norm for small residual norms. This is useful to conceal\n\ 2105ca161407SBarry Smith meaningless digits that may be different on different machines.\n",p); 210676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 2107ca161407SBarry Smith if (snes->type == SNES_NONLINEAR_EQUATIONS) { 210876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2109ca161407SBarry Smith " Options for solving systems of nonlinear equations only:\n"); 211076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 211176be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 211276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 211376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p); 211476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 211576be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2116ca161407SBarry Smith " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 211776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2118ca161407SBarry Smith " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 211976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2120ca161407SBarry Smith " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 212176be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2122ca161407SBarry Smith " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 212376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2124ca161407SBarry Smith " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 212576be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2126ca161407SBarry Smith " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 212776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2128ca161407SBarry Smith " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 2129ca161407SBarry Smith } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 213076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2131ca161407SBarry Smith " Options for solving unconstrained minimization problems only:\n"); 213276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 213376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 213476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 2135ca161407SBarry Smith } 213676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 213776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm,"a particular method\n"); 2138ca161407SBarry Smith if (snes->printhelp) { 2139ca161407SBarry Smith ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr); 2140ca161407SBarry Smith } 2141ca161407SBarry Smith PetscFunctionReturn(0); 2142ca161407SBarry Smith } 21433c7409f5SSatish Balay 2144acb85ae6SSatish Balay /*MC 2145b2002411SLois Curfman McInnes SNESRegister - Adds a method to the nonlinear solver package. 21469b94acceSBarry Smith 2147b2002411SLois Curfman McInnes Synopsis: 2148b2002411SLois Curfman McInnes SNESRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES)) 21499b94acceSBarry Smith 21508d76a1e5SLois Curfman McInnes Not collective 21518d76a1e5SLois Curfman McInnes 2152b2002411SLois Curfman McInnes Input Parameters: 2153c7afd0dbSLois Curfman McInnes + name_solver - name of a new user-defined solver 2154b2002411SLois Curfman McInnes . path - path (either absolute or relative) the library containing this solver 2155b2002411SLois Curfman McInnes . name_create - name of routine to create method context 2156c7afd0dbSLois Curfman McInnes - routine_create - routine to create method context 21579b94acceSBarry Smith 2158b2002411SLois Curfman McInnes Notes: 2159b2002411SLois Curfman McInnes SNESRegister() may be called multiple times to add several user-defined solvers. 21603a40ed3dSBarry Smith 2161b2002411SLois Curfman McInnes If dynamic libraries are used, then the fourth input argument (routine_create) 2162b2002411SLois Curfman McInnes is ignored. 2163b2002411SLois Curfman McInnes 2164b2002411SLois Curfman McInnes Sample usage: 216569225d78SLois Curfman McInnes .vb 2166b2002411SLois Curfman McInnes SNESRegister("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a, 2167b2002411SLois Curfman McInnes "MySolverCreate",MySolverCreate); 216869225d78SLois Curfman McInnes .ve 2169b2002411SLois Curfman McInnes 2170b2002411SLois Curfman McInnes Then, your solver can be chosen with the procedural interface via 2171b2002411SLois Curfman McInnes $ SNESSetType(snes,"my_solver") 2172b2002411SLois Curfman McInnes or at runtime via the option 2173b2002411SLois Curfman McInnes $ -snes_type my_solver 2174b2002411SLois Curfman McInnes 2175b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register 2176b2002411SLois Curfman McInnes 2177b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy() 2178b2002411SLois Curfman McInnes M*/ 2179b2002411SLois Curfman McInnes 2180b2002411SLois Curfman McInnes #undef __FUNC__ 2181b2002411SLois Curfman McInnes #define __FUNC__ "SNESRegister_Private" 2182b2002411SLois Curfman McInnes int SNESRegister_Private(char *sname,char *path,char *name,int (*function)(SNES)) 2183b2002411SLois Curfman McInnes { 2184b2002411SLois Curfman McInnes char fullname[256]; 2185b2002411SLois Curfman McInnes int ierr; 2186b2002411SLois Curfman McInnes 2187b2002411SLois Curfman McInnes PetscFunctionBegin; 2188b2002411SLois Curfman McInnes PetscStrcpy(fullname,path); PetscStrcat(fullname,":");PetscStrcat(fullname,name); 2189488ecbafSBarry Smith ierr = FListAdd_Private(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr); 2190b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2191b2002411SLois Curfman McInnes } 2192