15cd90555SBarry Smith 2a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 3*3e081fefSLois Curfman McInnes static char vcid[] = "$Id: snes.c,v 1.153 1998/04/24 22:11:25 curfman Exp curfman $"; 49b94acceSBarry Smith #endif 59b94acceSBarry Smith 670f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 7f5eb4b81SSatish Balay #include "src/sys/nreg.h" 89b94acceSBarry Smith #include "pinclude/pviewer.h" 99b94acceSBarry Smith #include <math.h> 109b94acceSBarry Smith 1184cb2905SBarry Smith int SNESRegisterAllCalled = 0; 1282bf6240SBarry Smith DLList SNESList = 0; 1382bf6240SBarry Smith 145615d1e5SSatish Balay #undef __FUNC__ 15d4bb536fSBarry Smith #define __FUNC__ "SNESView" 169b94acceSBarry Smith /*@ 179b94acceSBarry Smith SNESView - Prints the SNES data structure. 189b94acceSBarry Smith 19fee21e36SBarry Smith Collective on SNES, unless Viewer is VIEWER_STDOUT_SELF 20fee21e36SBarry Smith 21c7afd0dbSLois Curfman McInnes Input Parameters: 22c7afd0dbSLois Curfman McInnes + SNES - the SNES context 23c7afd0dbSLois Curfman McInnes - viewer - visualization context 24c7afd0dbSLois Curfman McInnes 259b94acceSBarry Smith Options Database Key: 26c8a8ba5cSLois Curfman McInnes . -snes_view - Calls SNESView() at end of SNESSolve() 279b94acceSBarry Smith 289b94acceSBarry Smith Notes: 299b94acceSBarry Smith The available visualization contexts include 30c7afd0dbSLois Curfman McInnes + VIEWER_STDOUT_SELF - standard output (default) 31c7afd0dbSLois Curfman McInnes - VIEWER_STDOUT_WORLD - synchronized standard 32c8a8ba5cSLois Curfman McInnes output where only the first processor opens 33c8a8ba5cSLois Curfman McInnes the file. All other processors send their 34c8a8ba5cSLois Curfman McInnes data to the first processor to print. 359b94acceSBarry Smith 36*3e081fefSLois Curfman McInnes The user can open an alternative visualization context with 371edb8f11SLois Curfman McInnes ViewerFileOpenASCII() - output to a specified file. 389b94acceSBarry Smith 399b94acceSBarry Smith .keywords: SNES, view 409b94acceSBarry Smith 41dbb450caSBarry Smith .seealso: ViewerFileOpenASCII() 429b94acceSBarry Smith @*/ 439b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 449b94acceSBarry Smith { 459b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 469b94acceSBarry Smith FILE *fd; 479b94acceSBarry Smith int ierr; 489b94acceSBarry Smith SLES sles; 499b94acceSBarry Smith char *method; 5019bcc07fSBarry Smith ViewerType vtype; 519b94acceSBarry Smith 523a40ed3dSBarry Smith PetscFunctionBegin; 5374679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5474679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5574679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5674679c65SBarry Smith 5719bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 5819bcc07fSBarry Smith if (vtype == ASCII_FILE_VIEWER || vtype == ASCII_FILES_VIEWER) { 5990ace30eSBarry Smith ierr = ViewerASCIIGetPointer(viewer,&fd); CHKERRQ(ierr); 6077c4ece6SBarry Smith PetscFPrintf(snes->comm,fd,"SNES Object:\n"); 6182bf6240SBarry Smith SNESGetType(snes,&method); 6277c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," method: %s\n",method); 63e1311b90SBarry Smith if (snes->view) (*snes->view)(snes,viewer); 6477c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 659b94acceSBarry Smith " maximum iterations=%d, maximum function evaluations=%d\n", 669b94acceSBarry Smith snes->max_its,snes->max_funcs); 6777c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 689b94acceSBarry Smith " tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 699b94acceSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol); 707a00f4a9SLois Curfman McInnes PetscFPrintf(snes->comm,fd, 717a00f4a9SLois Curfman McInnes " total number of linear solver iterations=%d\n",snes->linear_its); 72ae3c334cSLois Curfman McInnes PetscFPrintf(snes->comm,fd, 7368632166SLois Curfman McInnes " total number of function evaluations=%d\n",snes->nfuncs); 749b94acceSBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) 7577c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," min function tolerance=%g\n",snes->fmin); 769b94acceSBarry Smith if (snes->ksp_ewconv) { 779b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 789b94acceSBarry Smith if (kctx) { 7977c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 809b94acceSBarry Smith " Eisenstat-Walker computation of KSP relative tolerance (version %d)\n", 819b94acceSBarry Smith kctx->version); 8277c4ece6SBarry Smith PetscFPrintf(snes->comm,fd, 839b94acceSBarry Smith " rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0, 849b94acceSBarry Smith kctx->rtol_max,kctx->threshold); 8577c4ece6SBarry Smith PetscFPrintf(snes->comm,fd," gamma=%g, alpha=%g, alpha2=%g\n", 869b94acceSBarry Smith kctx->gamma,kctx->alpha,kctx->alpha2); 879b94acceSBarry Smith } 889b94acceSBarry Smith } 89c30f285eSLois Curfman McInnes } else if (vtype == STRING_VIEWER) { 9082bf6240SBarry Smith SNESGetType(snes,&method); 91c30f285eSLois Curfman McInnes ViewerStringSPrintf(viewer," %-3.3s",method); 9219bcc07fSBarry Smith } 939b94acceSBarry Smith SNESGetSLES(snes,&sles); 949b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 953a40ed3dSBarry Smith PetscFunctionReturn(0); 969b94acceSBarry Smith } 979b94acceSBarry Smith 98639f9d9dSBarry Smith /* 99639f9d9dSBarry Smith We retain a list of functions that also take SNES command 100639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 101639f9d9dSBarry Smith */ 102639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 103639f9d9dSBarry Smith static int numberofsetfromoptions; 104639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 105639f9d9dSBarry Smith 1065615d1e5SSatish Balay #undef __FUNC__ 107d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker" 108639f9d9dSBarry Smith /*@ 109639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 110639f9d9dSBarry Smith 111c7afd0dbSLois Curfman McInnes Not Collective 112c7afd0dbSLois Curfman McInnes 113639f9d9dSBarry Smith Input Parameter: 114639f9d9dSBarry Smith . snescheck - function that checks for options 115639f9d9dSBarry Smith 116639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 117639f9d9dSBarry Smith @*/ 118639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 119639f9d9dSBarry Smith { 1203a40ed3dSBarry Smith PetscFunctionBegin; 121639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 122a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed"); 123639f9d9dSBarry Smith } 124639f9d9dSBarry Smith 125639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 1263a40ed3dSBarry Smith PetscFunctionReturn(0); 127639f9d9dSBarry Smith } 128639f9d9dSBarry Smith 1295615d1e5SSatish Balay #undef __FUNC__ 1305615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1319b94acceSBarry Smith /*@ 132682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1339b94acceSBarry Smith 134c7afd0dbSLois Curfman McInnes Collective on SNES 135c7afd0dbSLois Curfman McInnes 1369b94acceSBarry Smith Input Parameter: 1379b94acceSBarry Smith . snes - the SNES context 1389b94acceSBarry Smith 13911ca99fdSLois Curfman McInnes Notes: 14011ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 14111ca99fdSLois Curfman McInnes the users manual. 14283e2fdc7SBarry Smith 1439b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 1449b94acceSBarry Smith 145a86d99e1SLois Curfman McInnes .seealso: SNESPrintHelp(), SNESSetOptionsPrefix() 1469b94acceSBarry Smith @*/ 1479b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 1489b94acceSBarry Smith { 14982bf6240SBarry Smith char method[256]; 1509b94acceSBarry Smith double tmp; 1519b94acceSBarry Smith SLES sles; 15281f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 1533c7409f5SSatish Balay int version = PETSC_DEFAULT; 1549b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 1559b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 1569b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 1579b94acceSBarry Smith double alpha = PETSC_DEFAULT; 1589b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 1599b94acceSBarry Smith double threshold = PETSC_DEFAULT; 1609b94acceSBarry Smith 1613a40ed3dSBarry Smith PetscFunctionBegin; 16281f57ec7SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 16381f57ec7SBarry Smith 16477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16582bf6240SBarry Smith if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp"); 166ca161407SBarry Smith 16782bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 16882bf6240SBarry Smith ierr = OptionsGetString(snes->prefix,"-snes_type",method,256,&flg); 169052efed2SBarry Smith if (flg) { 17082bf6240SBarry Smith ierr = SNESSetType(snes,(SNESType) method); CHKERRQ(ierr); 1715334005bSBarry Smith } 1723c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 173d64ed03dSBarry Smith if (flg) { 174d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 175d64ed03dSBarry Smith } 1763c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 177d64ed03dSBarry Smith if (flg) { 178d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 179d64ed03dSBarry Smith } 1803c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 181d64ed03dSBarry Smith if (flg) { 182d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 183d64ed03dSBarry Smith } 1843c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 1853c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 186d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 187d64ed03dSBarry Smith if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); } 188d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 189d64ed03dSBarry Smith if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp); CHKERRQ(ierr);} 1903c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 1913c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 192b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 193b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 194b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 195b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 196b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 197b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 198b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 19993c39befSBarry Smith 20093c39befSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr); 20193c39befSBarry Smith if (flg) {snes->converged = 0;} 20293c39befSBarry Smith 2039b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 2049b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 2055bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 2065cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 2073c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 208639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 2093c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 210d31a3109SLois Curfman McInnes if (flg) { 211639f9d9dSBarry Smith ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr); 212d31a3109SLois Curfman McInnes } 21381f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2143c7409f5SSatish Balay if (flg) { 21517699dbbSLois Curfman McInnes int rank = 0; 216d7e8b826SBarry Smith DrawLG lg; 21717699dbbSLois Curfman McInnes MPI_Initialized(&rank); 21817699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 21917699dbbSLois Curfman McInnes if (!rank) { 22081f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2219b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 222f8c826e1SBarry Smith snes->xmonitor = lg; 2239b94acceSBarry Smith PLogObjectParent(snes,lg); 2249b94acceSBarry Smith } 2259b94acceSBarry Smith } 2263c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2273c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2289b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2299b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 230981c4779SBarry Smith PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 231981c4779SBarry Smith } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 23231615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 23331615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 234d64ed03dSBarry Smith PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2359b94acceSBarry Smith } 236639f9d9dSBarry Smith 237639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 238639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 239639f9d9dSBarry Smith } 240639f9d9dSBarry Smith 2419b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 2429b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 24382bf6240SBarry Smith /* 24482bf6240SBarry Smith Since the private setfromoptions requires the type to all ready have 24582bf6240SBarry Smith been set we make sure a type is set by this time 24682bf6240SBarry Smith */ 24782bf6240SBarry Smith if (!snes->type_name) { 24882bf6240SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 24982bf6240SBarry Smith ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 25082bf6240SBarry Smith } else { 25182bf6240SBarry Smith ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 25282bf6240SBarry Smith } 25382bf6240SBarry Smith } 25482bf6240SBarry Smith 25582bf6240SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 25682bf6240SBarry Smith if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);} 25782bf6240SBarry Smith 2583a40ed3dSBarry Smith if (snes->setfromoptions) { 2593a40ed3dSBarry Smith ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 2603a40ed3dSBarry Smith } 2613a40ed3dSBarry Smith PetscFunctionReturn(0); 2629b94acceSBarry Smith } 2639b94acceSBarry Smith 264a847f771SSatish Balay 2655615d1e5SSatish Balay #undef __FUNC__ 266d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext" 2679b94acceSBarry Smith /*@ 2689b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 2699b94acceSBarry Smith the nonlinear solvers. 2709b94acceSBarry Smith 271fee21e36SBarry Smith Collective on SNES 272fee21e36SBarry Smith 273c7afd0dbSLois Curfman McInnes Input Parameters: 274c7afd0dbSLois Curfman McInnes + snes - the SNES context 275c7afd0dbSLois Curfman McInnes - usrP - optional user context 276c7afd0dbSLois Curfman McInnes 2779b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 2789b94acceSBarry Smith 2799b94acceSBarry Smith .seealso: SNESGetApplicationContext() 2809b94acceSBarry Smith @*/ 2819b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 2829b94acceSBarry Smith { 2833a40ed3dSBarry Smith PetscFunctionBegin; 28477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2859b94acceSBarry Smith snes->user = usrP; 2863a40ed3dSBarry Smith PetscFunctionReturn(0); 2879b94acceSBarry Smith } 28874679c65SBarry Smith 2895615d1e5SSatish Balay #undef __FUNC__ 290d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext" 2919b94acceSBarry Smith /*@C 2929b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 2939b94acceSBarry Smith nonlinear solvers. 2949b94acceSBarry Smith 295c7afd0dbSLois Curfman McInnes Not Collective 296c7afd0dbSLois Curfman McInnes 2979b94acceSBarry Smith Input Parameter: 2989b94acceSBarry Smith . snes - SNES context 2999b94acceSBarry Smith 3009b94acceSBarry Smith Output Parameter: 3019b94acceSBarry Smith . usrP - user context 3029b94acceSBarry Smith 3039b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3049b94acceSBarry Smith 3059b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3069b94acceSBarry Smith @*/ 3079b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3089b94acceSBarry Smith { 3093a40ed3dSBarry Smith PetscFunctionBegin; 31077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3119b94acceSBarry Smith *usrP = snes->user; 3123a40ed3dSBarry Smith PetscFunctionReturn(0); 3139b94acceSBarry Smith } 31474679c65SBarry Smith 3155615d1e5SSatish Balay #undef __FUNC__ 316d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber" 3179b94acceSBarry Smith /*@ 3189b94acceSBarry Smith SNESGetIterationNumber - Gets the current iteration number of the 3199b94acceSBarry Smith nonlinear solver. 3209b94acceSBarry Smith 321c7afd0dbSLois Curfman McInnes Not Collective 322c7afd0dbSLois Curfman McInnes 3239b94acceSBarry Smith Input Parameter: 3249b94acceSBarry Smith . snes - SNES context 3259b94acceSBarry Smith 3269b94acceSBarry Smith Output Parameter: 3279b94acceSBarry Smith . iter - iteration number 3289b94acceSBarry Smith 3299b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 3309b94acceSBarry Smith @*/ 3319b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 3329b94acceSBarry Smith { 3333a40ed3dSBarry Smith PetscFunctionBegin; 33477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 33574679c65SBarry Smith PetscValidIntPointer(iter); 3369b94acceSBarry Smith *iter = snes->iter; 3373a40ed3dSBarry Smith PetscFunctionReturn(0); 3389b94acceSBarry Smith } 33974679c65SBarry Smith 3405615d1e5SSatish Balay #undef __FUNC__ 3415615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 3429b94acceSBarry Smith /*@ 3439b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 3449b94acceSBarry Smith with SNESSSetFunction(). 3459b94acceSBarry Smith 346c7afd0dbSLois Curfman McInnes Collective on SNES 347c7afd0dbSLois Curfman McInnes 3489b94acceSBarry Smith Input Parameter: 3499b94acceSBarry Smith . snes - SNES context 3509b94acceSBarry Smith 3519b94acceSBarry Smith Output Parameter: 3529b94acceSBarry Smith . fnorm - 2-norm of function 3539b94acceSBarry Smith 3549b94acceSBarry Smith Note: 3559b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 356a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 357a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 3589b94acceSBarry Smith 3599b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 360a86d99e1SLois Curfman McInnes 361a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction() 3629b94acceSBarry Smith @*/ 3639b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 3649b94acceSBarry Smith { 3653a40ed3dSBarry Smith PetscFunctionBegin; 36677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 36774679c65SBarry Smith PetscValidScalarPointer(fnorm); 36874679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 369d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 37074679c65SBarry Smith } 3719b94acceSBarry Smith *fnorm = snes->norm; 3723a40ed3dSBarry Smith PetscFunctionReturn(0); 3739b94acceSBarry Smith } 37474679c65SBarry Smith 3755615d1e5SSatish Balay #undef __FUNC__ 3765615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 3779b94acceSBarry Smith /*@ 3789b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 3799b94acceSBarry Smith with SNESSSetGradient(). 3809b94acceSBarry Smith 381c7afd0dbSLois Curfman McInnes Collective on SNES 382c7afd0dbSLois Curfman McInnes 3839b94acceSBarry Smith Input Parameter: 3849b94acceSBarry Smith . snes - SNES context 3859b94acceSBarry Smith 3869b94acceSBarry Smith Output Parameter: 3879b94acceSBarry Smith . fnorm - 2-norm of gradient 3889b94acceSBarry Smith 3899b94acceSBarry Smith Note: 3909b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 391a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 392a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 3939b94acceSBarry Smith 3949b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 395a86d99e1SLois Curfman McInnes 396a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 3979b94acceSBarry Smith @*/ 3989b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 3999b94acceSBarry Smith { 4003a40ed3dSBarry Smith PetscFunctionBegin; 40177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 40274679c65SBarry Smith PetscValidScalarPointer(gnorm); 40374679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 404d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 40574679c65SBarry Smith } 4069b94acceSBarry Smith *gnorm = snes->norm; 4073a40ed3dSBarry Smith PetscFunctionReturn(0); 4089b94acceSBarry Smith } 40974679c65SBarry Smith 4105615d1e5SSatish Balay #undef __FUNC__ 411d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 4129b94acceSBarry Smith /*@ 4139b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4149b94acceSBarry Smith attempted by the nonlinear solver. 4159b94acceSBarry Smith 416c7afd0dbSLois Curfman McInnes Not Collective 417c7afd0dbSLois Curfman McInnes 4189b94acceSBarry Smith Input Parameter: 4199b94acceSBarry Smith . snes - SNES context 4209b94acceSBarry Smith 4219b94acceSBarry Smith Output Parameter: 4229b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 4239b94acceSBarry Smith 424c96a6f78SLois Curfman McInnes Notes: 425c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 426c96a6f78SLois Curfman McInnes 4279b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 4289b94acceSBarry Smith @*/ 4299b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 4309b94acceSBarry Smith { 4313a40ed3dSBarry Smith PetscFunctionBegin; 43277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 43374679c65SBarry Smith PetscValidIntPointer(nfails); 4349b94acceSBarry Smith *nfails = snes->nfailures; 4353a40ed3dSBarry Smith PetscFunctionReturn(0); 4369b94acceSBarry Smith } 437a847f771SSatish Balay 4385615d1e5SSatish Balay #undef __FUNC__ 439d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" 440c96a6f78SLois Curfman McInnes /*@ 441c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 442c96a6f78SLois Curfman McInnes used by the nonlinear solver. 443c96a6f78SLois Curfman McInnes 444c7afd0dbSLois Curfman McInnes Not Collective 445c7afd0dbSLois Curfman McInnes 446c96a6f78SLois Curfman McInnes Input Parameter: 447c96a6f78SLois Curfman McInnes . snes - SNES context 448c96a6f78SLois Curfman McInnes 449c96a6f78SLois Curfman McInnes Output Parameter: 450c96a6f78SLois Curfman McInnes . lits - number of linear iterations 451c96a6f78SLois Curfman McInnes 452c96a6f78SLois Curfman McInnes Notes: 453c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 454c96a6f78SLois Curfman McInnes 455c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 456c96a6f78SLois Curfman McInnes @*/ 45786bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 458c96a6f78SLois Curfman McInnes { 4593a40ed3dSBarry Smith PetscFunctionBegin; 460c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 461c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 462c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 4633a40ed3dSBarry Smith PetscFunctionReturn(0); 464c96a6f78SLois Curfman McInnes } 465c96a6f78SLois Curfman McInnes 466c96a6f78SLois Curfman McInnes #undef __FUNC__ 467d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES" 4689b94acceSBarry Smith /*@C 4699b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 4709b94acceSBarry Smith 471c7afd0dbSLois Curfman McInnes Not Collective, but if SNES object is parallel, then SLES object is parallel 472c7afd0dbSLois Curfman McInnes 4739b94acceSBarry Smith Input Parameter: 4749b94acceSBarry Smith . snes - the SNES context 4759b94acceSBarry Smith 4769b94acceSBarry Smith Output Parameter: 4779b94acceSBarry Smith . sles - the SLES context 4789b94acceSBarry Smith 4799b94acceSBarry Smith Notes: 4809b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 4819b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 4829b94acceSBarry Smith KSP and PC contexts as well. 4839b94acceSBarry Smith 4849b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 4859b94acceSBarry Smith 4869b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 4879b94acceSBarry Smith @*/ 4889b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 4899b94acceSBarry Smith { 4903a40ed3dSBarry Smith PetscFunctionBegin; 49177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 4929b94acceSBarry Smith *sles = snes->sles; 4933a40ed3dSBarry Smith PetscFunctionReturn(0); 4949b94acceSBarry Smith } 49582bf6240SBarry Smith 4969b94acceSBarry Smith /* -----------------------------------------------------------*/ 4975615d1e5SSatish Balay #undef __FUNC__ 4985615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 4999b94acceSBarry Smith /*@C 5009b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 5019b94acceSBarry Smith 502c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 503c7afd0dbSLois Curfman McInnes 504c7afd0dbSLois Curfman McInnes Input Parameters: 505c7afd0dbSLois Curfman McInnes + comm - MPI communicator 506c7afd0dbSLois Curfman McInnes - type - type of method, either 507c7afd0dbSLois Curfman McInnes SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations) 508c7afd0dbSLois Curfman McInnes or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization) 5099b94acceSBarry Smith 5109b94acceSBarry Smith Output Parameter: 5119b94acceSBarry Smith . outsnes - the new SNES context 5129b94acceSBarry Smith 513c7afd0dbSLois Curfman McInnes Options Database Keys: 514c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 515c7afd0dbSLois Curfman McInnes and no preconditioning matrix 516c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 517c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 518c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 519c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 520c1f60f51SBarry Smith 5219b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 5229b94acceSBarry Smith 52363a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 5249b94acceSBarry Smith @*/ 5254b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 5269b94acceSBarry Smith { 5279b94acceSBarry Smith int ierr; 5289b94acceSBarry Smith SNES snes; 5299b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 53037fcc0dbSBarry Smith 5313a40ed3dSBarry Smith PetscFunctionBegin; 5329b94acceSBarry Smith *outsnes = 0; 533d64ed03dSBarry Smith if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){ 534d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 535d64ed03dSBarry Smith } 536f830108cSBarry Smith PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,comm,SNESDestroy,SNESView); 5379b94acceSBarry Smith PLogObjectCreate(snes); 5389b94acceSBarry Smith snes->max_its = 50; 5399750a799SBarry Smith snes->max_funcs = 10000; 5409b94acceSBarry Smith snes->norm = 0.0; 5415a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 542b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 543b18e04deSLois Curfman McInnes snes->ttol = 0.0; 5449b94acceSBarry Smith snes->atol = 1.e-10; 545d64ed03dSBarry Smith } else { 546b4874afaSBarry Smith snes->rtol = 1.e-8; 547b4874afaSBarry Smith snes->ttol = 0.0; 548b4874afaSBarry Smith snes->atol = 1.e-50; 549b4874afaSBarry Smith } 5509b94acceSBarry Smith snes->xtol = 1.e-8; 551b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 5529b94acceSBarry Smith snes->nfuncs = 0; 5539b94acceSBarry Smith snes->nfailures = 0; 5547a00f4a9SLois Curfman McInnes snes->linear_its = 0; 555639f9d9dSBarry Smith snes->numbermonitors = 0; 5569b94acceSBarry Smith snes->data = 0; 5579b94acceSBarry Smith snes->view = 0; 5589b94acceSBarry Smith snes->computeumfunction = 0; 5599b94acceSBarry Smith snes->umfunP = 0; 5609b94acceSBarry Smith snes->fc = 0; 5619b94acceSBarry Smith snes->deltatol = 1.e-12; 5629b94acceSBarry Smith snes->fmin = -1.e30; 5639b94acceSBarry Smith snes->method_class = type; 5649b94acceSBarry Smith snes->set_method_called = 0; 56582bf6240SBarry Smith snes->setupcalled = 0; 5669b94acceSBarry Smith snes->ksp_ewconv = 0; 5676f24a144SLois Curfman McInnes snes->xmonitor = 0; 5686f24a144SLois Curfman McInnes snes->vwork = 0; 5696f24a144SLois Curfman McInnes snes->nwork = 0; 570c9005455SLois Curfman McInnes snes->conv_hist_size = 0; 571c9005455SLois Curfman McInnes snes->conv_act_size = 0; 572c9005455SLois Curfman McInnes snes->conv_hist = 0; 5739b94acceSBarry Smith 5749b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 5750452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 576eed86810SBarry Smith PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 5779b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 5789b94acceSBarry Smith kctx->version = 2; 5799b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 5809b94acceSBarry Smith this was too large for some test cases */ 5819b94acceSBarry Smith kctx->rtol_last = 0; 5829b94acceSBarry Smith kctx->rtol_max = .9; 5839b94acceSBarry Smith kctx->gamma = 1.0; 5849b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 5859b94acceSBarry Smith kctx->alpha = kctx->alpha2; 5869b94acceSBarry Smith kctx->threshold = .1; 5879b94acceSBarry Smith kctx->lresid_last = 0; 5889b94acceSBarry Smith kctx->norm_last = 0; 5899b94acceSBarry Smith 5909b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 5919b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 5925334005bSBarry Smith 5939b94acceSBarry Smith *outsnes = snes; 5943a40ed3dSBarry Smith PetscFunctionReturn(0); 5959b94acceSBarry Smith } 5969b94acceSBarry Smith 5979b94acceSBarry Smith /* --------------------------------------------------------------- */ 5985615d1e5SSatish Balay #undef __FUNC__ 599d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction" 6009b94acceSBarry Smith /*@C 6019b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 6029b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 6039b94acceSBarry Smith equations. 6049b94acceSBarry Smith 605fee21e36SBarry Smith Collective on SNES 606fee21e36SBarry Smith 607c7afd0dbSLois Curfman McInnes Input Parameters: 608c7afd0dbSLois Curfman McInnes + snes - the SNES context 609c7afd0dbSLois Curfman McInnes . func - function evaluation routine 610c7afd0dbSLois Curfman McInnes . r - vector to store function value 611c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 612c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 6139b94acceSBarry Smith 614c7afd0dbSLois Curfman McInnes Calling sequence of func: 615f40d9810SLois Curfman McInnes .vb 616f40d9810SLois Curfman McInnes func (SNES snes,Vec x,Vec f,void *ctx); 617f40d9810SLois Curfman McInnes .ve 618c7afd0dbSLois Curfman McInnes 619c7afd0dbSLois Curfman McInnes + x - input vector 620313e4042SLois Curfman McInnes . f - function vector 621c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 6229b94acceSBarry Smith 6239b94acceSBarry Smith Notes: 6249b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 6259b94acceSBarry Smith $ f'(x) x = -f(x), 626c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 6279b94acceSBarry Smith 6289b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6299b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6309b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 6319b94acceSBarry Smith 6329b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 6339b94acceSBarry Smith 634a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 6359b94acceSBarry Smith @*/ 6365334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 6379b94acceSBarry Smith { 6383a40ed3dSBarry Smith PetscFunctionBegin; 63977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 640a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 641a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 642a8c6a408SBarry Smith } 6439b94acceSBarry Smith snes->computefunction = func; 6449b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 6459b94acceSBarry Smith snes->funP = ctx; 6463a40ed3dSBarry Smith PetscFunctionReturn(0); 6479b94acceSBarry Smith } 6489b94acceSBarry Smith 6495615d1e5SSatish Balay #undef __FUNC__ 6505615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 6519b94acceSBarry Smith /*@ 6529b94acceSBarry Smith SNESComputeFunction - Computes the function that has been set with 6539b94acceSBarry Smith SNESSetFunction(). 6549b94acceSBarry Smith 655c7afd0dbSLois Curfman McInnes Collective on SNES 656c7afd0dbSLois Curfman McInnes 6579b94acceSBarry Smith Input Parameters: 658c7afd0dbSLois Curfman McInnes + snes - the SNES context 659c7afd0dbSLois Curfman McInnes - x - input vector 6609b94acceSBarry Smith 6619b94acceSBarry Smith Output Parameter: 6623638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 6639b94acceSBarry Smith 6641bffabb2SLois Curfman McInnes Notes: 6659b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 6669b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 6679b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 6689b94acceSBarry Smith 6699b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 6709b94acceSBarry Smith 671a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 6729b94acceSBarry Smith @*/ 6739b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 6749b94acceSBarry Smith { 6759b94acceSBarry Smith int ierr; 6769b94acceSBarry Smith 6773a40ed3dSBarry Smith PetscFunctionBegin; 678d4bb536fSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 679a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 680d4bb536fSBarry Smith } 6819b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 682d64ed03dSBarry Smith PetscStackPush("SNES user function"); 6839b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 684d64ed03dSBarry Smith PetscStackPop; 685ae3c334cSLois Curfman McInnes snes->nfuncs++; 6869b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 6873a40ed3dSBarry Smith PetscFunctionReturn(0); 6889b94acceSBarry Smith } 6899b94acceSBarry Smith 6905615d1e5SSatish Balay #undef __FUNC__ 691d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" 6929b94acceSBarry Smith /*@C 6939b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 6949b94acceSBarry Smith unconstrained minimization. 6959b94acceSBarry Smith 696fee21e36SBarry Smith Collective on SNES 697fee21e36SBarry Smith 698c7afd0dbSLois Curfman McInnes Input Parameters: 699c7afd0dbSLois Curfman McInnes + snes - the SNES context 700c7afd0dbSLois Curfman McInnes . func - function evaluation routine 701c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 702c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 7039b94acceSBarry Smith 704c7afd0dbSLois Curfman McInnes Calling sequence of func: 705f40d9810SLois Curfman McInnes .vb 706f40d9810SLois Curfman McInnes func (SNES snes,Vec x,double *f,void *ctx); 707f40d9810SLois Curfman McInnes .ve 708c7afd0dbSLois Curfman McInnes 709c7afd0dbSLois Curfman McInnes + x - input vector 7109b94acceSBarry Smith . f - function 711c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined function context 7129b94acceSBarry Smith 7139b94acceSBarry Smith Notes: 7149b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 7159b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 7169b94acceSBarry Smith SNESSetFunction(). 7179b94acceSBarry Smith 7189b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 7199b94acceSBarry Smith 720a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 721a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 7229b94acceSBarry Smith @*/ 7239b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 7249b94acceSBarry Smith void *ctx) 7259b94acceSBarry Smith { 7263a40ed3dSBarry Smith PetscFunctionBegin; 72777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 728a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 729a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 730a8c6a408SBarry Smith } 7319b94acceSBarry Smith snes->computeumfunction = func; 7329b94acceSBarry Smith snes->umfunP = ctx; 7333a40ed3dSBarry Smith PetscFunctionReturn(0); 7349b94acceSBarry Smith } 7359b94acceSBarry Smith 7365615d1e5SSatish Balay #undef __FUNC__ 7375615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 7389b94acceSBarry Smith /*@ 7399b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 7409b94acceSBarry Smith set with SNESSetMinimizationFunction(). 7419b94acceSBarry Smith 742c7afd0dbSLois Curfman McInnes Collective on SNES 743c7afd0dbSLois Curfman McInnes 7449b94acceSBarry Smith Input Parameters: 745c7afd0dbSLois Curfman McInnes + snes - the SNES context 746c7afd0dbSLois Curfman McInnes - x - input vector 7479b94acceSBarry Smith 7489b94acceSBarry Smith Output Parameter: 7499b94acceSBarry Smith . y - function value 7509b94acceSBarry Smith 7519b94acceSBarry Smith Notes: 7529b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 7539b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 7549b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 755a86d99e1SLois Curfman McInnes 756a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 757a86d99e1SLois Curfman McInnes 758a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 759a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 7609b94acceSBarry Smith @*/ 7619b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 7629b94acceSBarry Smith { 7639b94acceSBarry Smith int ierr; 7643a40ed3dSBarry Smith 7653a40ed3dSBarry Smith PetscFunctionBegin; 766a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 767a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 768a8c6a408SBarry Smith } 7699b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 770d64ed03dSBarry Smith PetscStackPush("SNES user minimzation function"); 7719b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 772d64ed03dSBarry Smith PetscStackPop; 773ae3c334cSLois Curfman McInnes snes->nfuncs++; 7749b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 7753a40ed3dSBarry Smith PetscFunctionReturn(0); 7769b94acceSBarry Smith } 7779b94acceSBarry Smith 7785615d1e5SSatish Balay #undef __FUNC__ 779d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient" 7809b94acceSBarry Smith /*@C 7819b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 7829b94acceSBarry Smith vector for use by the SNES routines. 7839b94acceSBarry Smith 784c7afd0dbSLois Curfman McInnes Collective on SNES 785c7afd0dbSLois Curfman McInnes 7869b94acceSBarry Smith Input Parameters: 787c7afd0dbSLois Curfman McInnes + snes - the SNES context 7889b94acceSBarry Smith . func - function evaluation routine 789044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 790044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 791c7afd0dbSLois Curfman McInnes - r - vector to store gradient value 792fee21e36SBarry Smith 7939b94acceSBarry Smith Calling sequence of func: 794f40d9810SLois Curfman McInnes .vb 795f40d9810SLois Curfman McInnes func (SNES, Vec x, Vec g, void *ctx); 796f40d9810SLois Curfman McInnes .ve 7979b94acceSBarry Smith 798c7afd0dbSLois Curfman McInnes + x - input vector 7999b94acceSBarry Smith . g - gradient vector 800c7afd0dbSLois Curfman McInnes - ctx - optional user-defined gradient context 8019b94acceSBarry Smith 8029b94acceSBarry Smith Notes: 8039b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 8049b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 8059b94acceSBarry Smith SNESSetFunction(). 8069b94acceSBarry Smith 8079b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 8089b94acceSBarry Smith 809a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 810a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 8119b94acceSBarry Smith @*/ 81274679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 8139b94acceSBarry Smith { 8143a40ed3dSBarry Smith PetscFunctionBegin; 81577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 816a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 817a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 818a8c6a408SBarry Smith } 8199b94acceSBarry Smith snes->computefunction = func; 8209b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 8219b94acceSBarry Smith snes->funP = ctx; 8223a40ed3dSBarry Smith PetscFunctionReturn(0); 8239b94acceSBarry Smith } 8249b94acceSBarry Smith 8255615d1e5SSatish Balay #undef __FUNC__ 8265615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 8279b94acceSBarry Smith /*@ 828a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 829a86d99e1SLois Curfman McInnes SNESSetGradient(). 8309b94acceSBarry Smith 831c7afd0dbSLois Curfman McInnes Collective on SNES 832c7afd0dbSLois Curfman McInnes 8339b94acceSBarry Smith Input Parameters: 834c7afd0dbSLois Curfman McInnes + snes - the SNES context 835c7afd0dbSLois Curfman McInnes - x - input vector 8369b94acceSBarry Smith 8379b94acceSBarry Smith Output Parameter: 8389b94acceSBarry Smith . y - gradient vector 8399b94acceSBarry Smith 8409b94acceSBarry Smith Notes: 8419b94acceSBarry Smith SNESComputeGradient() is valid only for 8429b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 8439b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 844a86d99e1SLois Curfman McInnes 845a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 846a86d99e1SLois Curfman McInnes 847a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 848a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 8499b94acceSBarry Smith @*/ 8509b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 8519b94acceSBarry Smith { 8529b94acceSBarry Smith int ierr; 8533a40ed3dSBarry Smith 8543a40ed3dSBarry Smith PetscFunctionBegin; 8553a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 856a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 8573a40ed3dSBarry Smith } 8589b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 859d64ed03dSBarry Smith PetscStackPush("SNES user gradient function"); 8609b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 861d64ed03dSBarry Smith PetscStackPop; 8629b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 8633a40ed3dSBarry Smith PetscFunctionReturn(0); 8649b94acceSBarry Smith } 8659b94acceSBarry Smith 8665615d1e5SSatish Balay #undef __FUNC__ 8675615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 86862fef451SLois Curfman McInnes /*@ 86962fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 87062fef451SLois Curfman McInnes set with SNESSetJacobian(). 87162fef451SLois Curfman McInnes 872c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 873c7afd0dbSLois Curfman McInnes 87462fef451SLois Curfman McInnes Input Parameters: 875c7afd0dbSLois Curfman McInnes + snes - the SNES context 876c7afd0dbSLois Curfman McInnes - x - input vector 87762fef451SLois Curfman McInnes 87862fef451SLois Curfman McInnes Output Parameters: 879c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 88062fef451SLois Curfman McInnes . B - optional preconditioning matrix 881c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 882fee21e36SBarry Smith 88362fef451SLois Curfman McInnes Notes: 88462fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 88562fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 88662fef451SLois Curfman McInnes 887dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 888dc5a77f8SLois Curfman McInnes flag parameter. 88962fef451SLois Curfman McInnes 89062fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 89162fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 892005c665bSBarry Smith methods is SNESComputeHessian(). 89362fef451SLois Curfman McInnes 89462fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 89562fef451SLois Curfman McInnes 89662fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 89762fef451SLois Curfman McInnes @*/ 8989b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 8999b94acceSBarry Smith { 9009b94acceSBarry Smith int ierr; 9013a40ed3dSBarry Smith 9023a40ed3dSBarry Smith PetscFunctionBegin; 9033a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 904a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 9053a40ed3dSBarry Smith } 9063a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 9079b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 908c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 909d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 9109b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 911d64ed03dSBarry Smith PetscStackPop; 9129b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 9136d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 91477c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 91577c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9163a40ed3dSBarry Smith PetscFunctionReturn(0); 9179b94acceSBarry Smith } 9189b94acceSBarry Smith 9195615d1e5SSatish Balay #undef __FUNC__ 9205615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 92162fef451SLois Curfman McInnes /*@ 92262fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 92362fef451SLois Curfman McInnes set with SNESSetHessian(). 92462fef451SLois Curfman McInnes 925c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 926c7afd0dbSLois Curfman McInnes 92762fef451SLois Curfman McInnes Input Parameters: 928c7afd0dbSLois Curfman McInnes + snes - the SNES context 929c7afd0dbSLois Curfman McInnes - x - input vector 93062fef451SLois Curfman McInnes 93162fef451SLois Curfman McInnes Output Parameters: 932c7afd0dbSLois Curfman McInnes + A - Hessian matrix 93362fef451SLois Curfman McInnes . B - optional preconditioning matrix 934c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 935fee21e36SBarry Smith 93662fef451SLois Curfman McInnes Notes: 93762fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 93862fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 93962fef451SLois Curfman McInnes 940dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 941dc5a77f8SLois Curfman McInnes flag parameter. 94262fef451SLois Curfman McInnes 94362fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 94462fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 94562fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 94662fef451SLois Curfman McInnes 94762fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 94862fef451SLois Curfman McInnes 949a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 950a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 95162fef451SLois Curfman McInnes @*/ 95262fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 9539b94acceSBarry Smith { 9549b94acceSBarry Smith int ierr; 9553a40ed3dSBarry Smith 9563a40ed3dSBarry Smith PetscFunctionBegin; 9573a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 958a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 9593a40ed3dSBarry Smith } 9603a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 96162fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 9620f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 963d64ed03dSBarry Smith PetscStackPush("SNES user Hessian function"); 96462fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 965d64ed03dSBarry Smith PetscStackPop; 96662fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 9670f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 96877c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 96977c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 9703a40ed3dSBarry Smith PetscFunctionReturn(0); 9719b94acceSBarry Smith } 9729b94acceSBarry Smith 9735615d1e5SSatish Balay #undef __FUNC__ 974d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian" 9759b94acceSBarry Smith /*@C 9769b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 977044dda88SLois Curfman McInnes location to store the matrix. 9789b94acceSBarry Smith 979c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 980c7afd0dbSLois Curfman McInnes 9819b94acceSBarry Smith Input Parameters: 982c7afd0dbSLois Curfman McInnes + snes - the SNES context 9839b94acceSBarry Smith . A - Jacobian matrix 9849b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 9859b94acceSBarry Smith . func - Jacobian evaluation routine 986c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 9872cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 9889b94acceSBarry Smith 9899b94acceSBarry Smith Calling sequence of func: 990f40d9810SLois Curfman McInnes .vb 991f40d9810SLois Curfman McInnes func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 992f40d9810SLois Curfman McInnes .ve 9939b94acceSBarry Smith 994c7afd0dbSLois Curfman McInnes + x - input vector 9959b94acceSBarry Smith . A - Jacobian matrix 9969b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 997ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 998ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 999c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 10009b94acceSBarry Smith 10019b94acceSBarry Smith Notes: 1002dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10032cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1004ac21db08SLois Curfman McInnes 1005ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 10069b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 10079b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10089b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10099b94acceSBarry Smith throughout the global iterations. 10109b94acceSBarry Smith 10119b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 10129b94acceSBarry Smith 1013ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 10149b94acceSBarry Smith @*/ 10159b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 10169b94acceSBarry Smith MatStructure*,void*),void *ctx) 10179b94acceSBarry Smith { 10183a40ed3dSBarry Smith PetscFunctionBegin; 101977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1020a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1021a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1022a8c6a408SBarry Smith } 10239b94acceSBarry Smith snes->computejacobian = func; 10249b94acceSBarry Smith snes->jacP = ctx; 10259b94acceSBarry Smith snes->jacobian = A; 10269b94acceSBarry Smith snes->jacobian_pre = B; 10273a40ed3dSBarry Smith PetscFunctionReturn(0); 10289b94acceSBarry Smith } 102962fef451SLois Curfman McInnes 10305615d1e5SSatish Balay #undef __FUNC__ 1031d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian" 1032b4fd4287SBarry Smith /*@ 1033b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1034b4fd4287SBarry Smith provided context for evaluating the Jacobian. 1035b4fd4287SBarry Smith 1036c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 1037c7afd0dbSLois Curfman McInnes 1038b4fd4287SBarry Smith Input Parameter: 1039b4fd4287SBarry Smith . snes - the nonlinear solver context 1040b4fd4287SBarry Smith 1041b4fd4287SBarry Smith Output Parameters: 1042c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 1043b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 1044c7afd0dbSLois Curfman McInnes - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1045fee21e36SBarry Smith 1046b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 1047b4fd4287SBarry Smith @*/ 1048b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1049b4fd4287SBarry Smith { 10503a40ed3dSBarry Smith PetscFunctionBegin; 10513a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1052a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 10533a40ed3dSBarry Smith } 1054b4fd4287SBarry Smith if (A) *A = snes->jacobian; 1055b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 1056b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 10573a40ed3dSBarry Smith PetscFunctionReturn(0); 1058b4fd4287SBarry Smith } 1059b4fd4287SBarry Smith 10605615d1e5SSatish Balay #undef __FUNC__ 1061d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian" 10629b94acceSBarry Smith /*@C 10639b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 1064044dda88SLois Curfman McInnes location to store the matrix. 10659b94acceSBarry Smith 1066c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1067c7afd0dbSLois Curfman McInnes 10689b94acceSBarry Smith Input Parameters: 1069c7afd0dbSLois Curfman McInnes + snes - the SNES context 10709b94acceSBarry Smith . A - Hessian matrix 10719b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 10729b94acceSBarry Smith . func - Jacobian evaluation routine 1073c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 1074313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 10759b94acceSBarry Smith 10769b94acceSBarry Smith Calling sequence of func: 1077f40d9810SLois Curfman McInnes .vb 1078f40d9810SLois Curfman McInnes func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1079f40d9810SLois Curfman McInnes .ve 10809b94acceSBarry Smith 1081c7afd0dbSLois Curfman McInnes + x - input vector 10829b94acceSBarry Smith . A - Hessian matrix 10839b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1084ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1085ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 1086c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Hessian context 10879b94acceSBarry Smith 10889b94acceSBarry Smith Notes: 1089dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 10902cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1091ac21db08SLois Curfman McInnes 10929b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 10939b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 10949b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 10959b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 10969b94acceSBarry Smith throughout the global iterations. 10979b94acceSBarry Smith 10989b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 10999b94acceSBarry Smith 1100ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 11019b94acceSBarry Smith @*/ 11029b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 11039b94acceSBarry Smith MatStructure*,void*),void *ctx) 11049b94acceSBarry Smith { 11053a40ed3dSBarry Smith PetscFunctionBegin; 110677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1107d4bb536fSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1108a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1109d4bb536fSBarry Smith } 11109b94acceSBarry Smith snes->computejacobian = func; 11119b94acceSBarry Smith snes->jacP = ctx; 11129b94acceSBarry Smith snes->jacobian = A; 11139b94acceSBarry Smith snes->jacobian_pre = B; 11143a40ed3dSBarry Smith PetscFunctionReturn(0); 11159b94acceSBarry Smith } 11169b94acceSBarry Smith 11175615d1e5SSatish Balay #undef __FUNC__ 1118d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian" 111962fef451SLois Curfman McInnes /*@ 112062fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 112162fef451SLois Curfman McInnes provided context for evaluating the Hessian. 112262fef451SLois Curfman McInnes 1123c7afd0dbSLois Curfman McInnes Not Collective, but Mat object is parallel if SNES object is parallel 1124c7afd0dbSLois Curfman McInnes 112562fef451SLois Curfman McInnes Input Parameter: 112662fef451SLois Curfman McInnes . snes - the nonlinear solver context 112762fef451SLois Curfman McInnes 112862fef451SLois Curfman McInnes Output Parameters: 1129c7afd0dbSLois Curfman McInnes + A - location to stash Hessian matrix (or PETSC_NULL) 113062fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 1131c7afd0dbSLois Curfman McInnes - ctx - location to stash Hessian ctx (or PETSC_NULL) 1132fee21e36SBarry Smith 113362fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 1134c7afd0dbSLois Curfman McInnes 1135c7afd0dbSLois Curfman McInnes .keywords: SNES, get, Hessian 113662fef451SLois Curfman McInnes @*/ 113762fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 113862fef451SLois Curfman McInnes { 11393a40ed3dSBarry Smith PetscFunctionBegin; 11403a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){ 1141a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 11423a40ed3dSBarry Smith } 114362fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 114462fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 114562fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 11463a40ed3dSBarry Smith PetscFunctionReturn(0); 114762fef451SLois Curfman McInnes } 114862fef451SLois Curfman McInnes 11499b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 11509b94acceSBarry Smith 11515615d1e5SSatish Balay #undef __FUNC__ 11525615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 11539b94acceSBarry Smith /*@ 11549b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1155272ac6f2SLois Curfman McInnes of a nonlinear solver. 11569b94acceSBarry Smith 1157fee21e36SBarry Smith Collective on SNES 1158fee21e36SBarry Smith 1159c7afd0dbSLois Curfman McInnes Input Parameters: 1160c7afd0dbSLois Curfman McInnes + snes - the SNES context 1161c7afd0dbSLois Curfman McInnes - x - the solution vector 1162c7afd0dbSLois Curfman McInnes 1163272ac6f2SLois Curfman McInnes Notes: 1164272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1165272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1166272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1167272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1168272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1169272ac6f2SLois Curfman McInnes 11709b94acceSBarry Smith .keywords: SNES, nonlinear, setup 11719b94acceSBarry Smith 11729b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 11739b94acceSBarry Smith @*/ 11748ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 11759b94acceSBarry Smith { 1176272ac6f2SLois Curfman McInnes int ierr, flg; 11773a40ed3dSBarry Smith 11783a40ed3dSBarry Smith PetscFunctionBegin; 117977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 118077c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 11818ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 11829b94acceSBarry Smith 1183c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1184c1f60f51SBarry Smith /* 1185c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1186dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1187c1f60f51SBarry Smith */ 1188c1f60f51SBarry Smith if (flg) { 1189c1f60f51SBarry Smith Mat J; 1190c1f60f51SBarry Smith ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1191c1f60f51SBarry Smith PLogObjectParent(snes,J); 1192c1f60f51SBarry Smith snes->mfshell = J; 1193c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1194c1f60f51SBarry Smith snes->jacobian = J; 1195c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1196d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1197c1f60f51SBarry Smith snes->jacobian = J; 1198c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1199d4bb536fSBarry Smith } else { 1200a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option"); 1201d4bb536fSBarry Smith } 1202c1f60f51SBarry Smith } 1203272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1204c1f60f51SBarry Smith /* 1205dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1206c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1207c1f60f51SBarry Smith */ 120831615d04SLois Curfman McInnes if (flg) { 1209272ac6f2SLois Curfman McInnes Mat J; 1210272ac6f2SLois Curfman McInnes ierr = SNESDefaultMatrixFreeMatCreate(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1211272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1212272ac6f2SLois Curfman McInnes snes->mfshell = J; 121383e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 121483e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 121594a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1216d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 121783e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 121894a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1219d4bb536fSBarry Smith } else { 1220a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option"); 1221d4bb536fSBarry Smith } 1222272ac6f2SLois Curfman McInnes } 12239b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1224a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1225a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1226a8c6a408SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first"); 1227a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 1228a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector"); 1229a8c6a408SBarry Smith } 1230a703fe33SLois Curfman McInnes 1231a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 123282bf6240SBarry Smith if (snes->ksp_ewconv && PetscStrcmp(snes->type_name,SNES_EQ_TR)) { 1233a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1234a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1235a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 1236a703fe33SLois Curfman McInnes ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private, 1237a703fe33SLois Curfman McInnes (void *)snes); CHKERRQ(ierr); 1238a703fe33SLois Curfman McInnes } 12399b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1240a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1241a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1242a8c6a408SBarry Smith if (!snes->computeumfunction) { 1243a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first"); 1244a8c6a408SBarry Smith } 1245a8c6a408SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian() first"); 1246d4bb536fSBarry Smith } else { 1247a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class"); 1248d4bb536fSBarry Smith } 1249a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 125082bf6240SBarry Smith snes->setupcalled = 1; 12513a40ed3dSBarry Smith PetscFunctionReturn(0); 12529b94acceSBarry Smith } 12539b94acceSBarry Smith 12545615d1e5SSatish Balay #undef __FUNC__ 1255d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy" 12569b94acceSBarry Smith /*@C 12579b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 12589b94acceSBarry Smith with SNESCreate(). 12599b94acceSBarry Smith 1260c7afd0dbSLois Curfman McInnes Collective on SNES 1261c7afd0dbSLois Curfman McInnes 12629b94acceSBarry Smith Input Parameter: 12639b94acceSBarry Smith . snes - the SNES context 12649b94acceSBarry Smith 12659b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 12669b94acceSBarry Smith 126763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 12689b94acceSBarry Smith @*/ 12699b94acceSBarry Smith int SNESDestroy(SNES snes) 12709b94acceSBarry Smith { 12719b94acceSBarry Smith int ierr; 12723a40ed3dSBarry Smith 12733a40ed3dSBarry Smith PetscFunctionBegin; 127477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12753a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1276d4bb536fSBarry Smith 1277e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes); CHKERRQ(ierr);} 12780452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 1279522c5e43SBarry Smith if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);} 12809b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 1281522c5e43SBarry Smith if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);} 1282522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 12839b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 12840452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 12853a40ed3dSBarry Smith PetscFunctionReturn(0); 12869b94acceSBarry Smith } 12879b94acceSBarry Smith 12889b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 12899b94acceSBarry Smith 12905615d1e5SSatish Balay #undef __FUNC__ 12915615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 12929b94acceSBarry Smith /*@ 1293d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 12949b94acceSBarry Smith 1295c7afd0dbSLois Curfman McInnes Collective on SNES 1296c7afd0dbSLois Curfman McInnes 12979b94acceSBarry Smith Input Parameters: 1298c7afd0dbSLois Curfman McInnes + snes - the SNES context 129933174efeSLois Curfman McInnes . atol - absolute convergence tolerance 130033174efeSLois Curfman McInnes . rtol - relative convergence tolerance 130133174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 130233174efeSLois Curfman McInnes of the change in the solution between steps 130333174efeSLois Curfman McInnes . maxit - maximum number of iterations 1304c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1305fee21e36SBarry Smith 130633174efeSLois Curfman McInnes Options Database Keys: 1307c7afd0dbSLois Curfman McInnes + -snes_atol <atol> - Sets atol 1308c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1309c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1310c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1311c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 13129b94acceSBarry Smith 1313d7a720efSLois Curfman McInnes Notes: 13149b94acceSBarry Smith The default maximum number of iterations is 50. 13159b94acceSBarry Smith The default maximum number of function evaluations is 1000. 13169b94acceSBarry Smith 131733174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 13189b94acceSBarry Smith 1319d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 13209b94acceSBarry Smith @*/ 1321d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 13229b94acceSBarry Smith { 13233a40ed3dSBarry Smith PetscFunctionBegin; 132477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1325d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1326d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1327d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1328d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1329d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 13303a40ed3dSBarry Smith PetscFunctionReturn(0); 13319b94acceSBarry Smith } 13329b94acceSBarry Smith 13335615d1e5SSatish Balay #undef __FUNC__ 13345615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 13359b94acceSBarry Smith /*@ 133633174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 133733174efeSLois Curfman McInnes 1338c7afd0dbSLois Curfman McInnes Not Collective 1339c7afd0dbSLois Curfman McInnes 134033174efeSLois Curfman McInnes Input Parameters: 1341c7afd0dbSLois Curfman McInnes + snes - the SNES context 134233174efeSLois Curfman McInnes . atol - absolute convergence tolerance 134333174efeSLois Curfman McInnes . rtol - relative convergence tolerance 134433174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 134533174efeSLois Curfman McInnes of the change in the solution between steps 134633174efeSLois Curfman McInnes . maxit - maximum number of iterations 1347c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1348fee21e36SBarry Smith 134933174efeSLois Curfman McInnes Notes: 135033174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 135133174efeSLois Curfman McInnes 135233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 135333174efeSLois Curfman McInnes 135433174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 135533174efeSLois Curfman McInnes @*/ 135633174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 135733174efeSLois Curfman McInnes { 13583a40ed3dSBarry Smith PetscFunctionBegin; 135933174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 136033174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 136133174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 136233174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 136333174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 136433174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 13653a40ed3dSBarry Smith PetscFunctionReturn(0); 136633174efeSLois Curfman McInnes } 136733174efeSLois Curfman McInnes 13685615d1e5SSatish Balay #undef __FUNC__ 13695615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 137033174efeSLois Curfman McInnes /*@ 13719b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 13729b94acceSBarry Smith 1373fee21e36SBarry Smith Collective on SNES 1374fee21e36SBarry Smith 1375c7afd0dbSLois Curfman McInnes Input Parameters: 1376c7afd0dbSLois Curfman McInnes + snes - the SNES context 1377c7afd0dbSLois Curfman McInnes - tol - tolerance 1378c7afd0dbSLois Curfman McInnes 13799b94acceSBarry Smith Options Database Key: 1380c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 13819b94acceSBarry Smith 13829b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 13839b94acceSBarry Smith 1384d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 13859b94acceSBarry Smith @*/ 13869b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 13879b94acceSBarry Smith { 13883a40ed3dSBarry Smith PetscFunctionBegin; 138977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13909b94acceSBarry Smith snes->deltatol = tol; 13913a40ed3dSBarry Smith PetscFunctionReturn(0); 13929b94acceSBarry Smith } 13939b94acceSBarry Smith 13945615d1e5SSatish Balay #undef __FUNC__ 13955615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 13969b94acceSBarry Smith /*@ 139777c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 13989b94acceSBarry Smith for unconstrained minimization solvers. 13999b94acceSBarry Smith 1400fee21e36SBarry Smith Collective on SNES 1401fee21e36SBarry Smith 1402c7afd0dbSLois Curfman McInnes Input Parameters: 1403c7afd0dbSLois Curfman McInnes + snes - the SNES context 1404c7afd0dbSLois Curfman McInnes - ftol - minimum function tolerance 1405c7afd0dbSLois Curfman McInnes 14069b94acceSBarry Smith Options Database Key: 1407c7afd0dbSLois Curfman McInnes . -snes_fmin <ftol> - Sets ftol 14089b94acceSBarry Smith 14099b94acceSBarry Smith Note: 141077c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 14119b94acceSBarry Smith methods only. 14129b94acceSBarry Smith 14139b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 14149b94acceSBarry Smith 1415d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 14169b94acceSBarry Smith @*/ 141777c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 14189b94acceSBarry Smith { 14193a40ed3dSBarry Smith PetscFunctionBegin; 142077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 14219b94acceSBarry Smith snes->fmin = ftol; 14223a40ed3dSBarry Smith PetscFunctionReturn(0); 14239b94acceSBarry Smith } 14249b94acceSBarry Smith 14259b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 14269b94acceSBarry Smith 14275615d1e5SSatish Balay #undef __FUNC__ 1428d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor" 14299b94acceSBarry Smith /*@C 14305cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 14319b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 14329b94acceSBarry Smith progress. 14339b94acceSBarry Smith 1434fee21e36SBarry Smith Collective on SNES 1435fee21e36SBarry Smith 1436c7afd0dbSLois Curfman McInnes Input Parameters: 1437c7afd0dbSLois Curfman McInnes + snes - the SNES context 1438c7afd0dbSLois Curfman McInnes . func - monitoring routine 1439c7afd0dbSLois Curfman McInnes - mctx - [optional] user-defined context for private data for the 1440c7afd0dbSLois Curfman McInnes monitor routine (may be PETSC_NULL) 14419b94acceSBarry Smith 1442c7afd0dbSLois Curfman McInnes Calling sequence of func: 1443f40d9810SLois Curfman McInnes .vb 1444f40d9810SLois Curfman McInnes int func(SNES snes,int its, Vec x,Vec f,double norm,void *mctx) 1445f40d9810SLois Curfman McInnes .ve 1446c7afd0dbSLois Curfman McInnes 1447c7afd0dbSLois Curfman McInnes + snes - the SNES context 1448c7afd0dbSLois Curfman McInnes . its - iteration number 1449c7afd0dbSLois Curfman McInnes . mctx - [optional] monitoring context 1450c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 1451c7afd0dbSLois Curfman McInnes for SNES_NONLINEAR_EQUATIONS methods 1452c7afd0dbSLois Curfman McInnes - norm - 2-norm gradient value (may be estimated) 1453c7afd0dbSLois Curfman McInnes for SNES_UNCONSTRAINED_MINIMIZATION methods 14549b94acceSBarry Smith 14559665c990SLois Curfman McInnes Options Database Keys: 1456c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1457c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1458c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1459c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1460c7afd0dbSLois Curfman McInnes been hardwired into a code by 1461c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1462c7afd0dbSLois Curfman McInnes does not cancel those set via 1463c7afd0dbSLois Curfman McInnes the options database. 14649665c990SLois Curfman McInnes 1465639f9d9dSBarry Smith Notes: 14666bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 14676bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 14686bc08f3fSLois Curfman McInnes order in which they were set. 1469639f9d9dSBarry Smith 14709b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 14719b94acceSBarry Smith 14725cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 14739b94acceSBarry Smith @*/ 147474679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 14759b94acceSBarry Smith { 14763a40ed3dSBarry Smith PetscFunctionBegin; 1477639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1478a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set"); 1479639f9d9dSBarry Smith } 1480639f9d9dSBarry Smith 1481639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1482639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 14833a40ed3dSBarry Smith PetscFunctionReturn(0); 14849b94acceSBarry Smith } 14859b94acceSBarry Smith 14865615d1e5SSatish Balay #undef __FUNC__ 14875cd90555SBarry Smith #define __FUNC__ "SNESClearMonitor" 14885cd90555SBarry Smith /*@C 14895cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 14905cd90555SBarry Smith 1491c7afd0dbSLois Curfman McInnes Collective on SNES 1492c7afd0dbSLois Curfman McInnes 14935cd90555SBarry Smith Input Parameters: 14945cd90555SBarry Smith . snes - the SNES context 14955cd90555SBarry Smith 14965cd90555SBarry Smith Options Database: 1497c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1498c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1499c7afd0dbSLois Curfman McInnes set via the options database 15005cd90555SBarry Smith 15015cd90555SBarry Smith Notes: 15025cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 15035cd90555SBarry Smith 15045cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 15055cd90555SBarry Smith 15065cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 15075cd90555SBarry Smith @*/ 15085cd90555SBarry Smith int SNESClearMonitor( SNES snes ) 15095cd90555SBarry Smith { 15105cd90555SBarry Smith PetscFunctionBegin; 15115cd90555SBarry Smith snes->numbermonitors = 0; 15125cd90555SBarry Smith PetscFunctionReturn(0); 15135cd90555SBarry Smith } 15145cd90555SBarry Smith 15155cd90555SBarry Smith #undef __FUNC__ 1516d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest" 15179b94acceSBarry Smith /*@C 15189b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 15199b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 15209b94acceSBarry Smith 1521fee21e36SBarry Smith Collective on SNES 1522fee21e36SBarry Smith 1523c7afd0dbSLois Curfman McInnes Input Parameters: 1524c7afd0dbSLois Curfman McInnes + snes - the SNES context 1525c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1526c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1527c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 15289b94acceSBarry Smith 1529c7afd0dbSLois Curfman McInnes Calling sequence of func: 1530f40d9810SLois Curfman McInnes .vb 1531f40d9810SLois Curfman McInnes int func (SNES snes,double xnorm,double gnorm,double f,void *cctx) 1532f40d9810SLois Curfman McInnes .ve 1533c7afd0dbSLois Curfman McInnes 1534c7afd0dbSLois Curfman McInnes + snes - the SNES context 1535c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1536c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 1537c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods) 1538c7afd0dbSLois Curfman McInnes . f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods) 1539c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods) 1540c7afd0dbSLois Curfman McInnes - f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods) 15419b94acceSBarry Smith 15429b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 15439b94acceSBarry Smith 154440191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 154540191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 15469b94acceSBarry Smith @*/ 154774679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 15489b94acceSBarry Smith { 15493a40ed3dSBarry Smith PetscFunctionBegin; 15509b94acceSBarry Smith (snes)->converged = func; 15519b94acceSBarry Smith (snes)->cnvP = cctx; 15523a40ed3dSBarry Smith PetscFunctionReturn(0); 15539b94acceSBarry Smith } 15549b94acceSBarry Smith 15555615d1e5SSatish Balay #undef __FUNC__ 1556d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" 1557c9005455SLois Curfman McInnes /*@ 1558c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1559c9005455SLois Curfman McInnes 1560fee21e36SBarry Smith Collective on SNES 1561fee21e36SBarry Smith 1562c7afd0dbSLois Curfman McInnes Input Parameters: 1563c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1564c7afd0dbSLois Curfman McInnes . a - array to hold history 1565c7afd0dbSLois Curfman McInnes - na - size of a 1566c7afd0dbSLois Curfman McInnes 1567c9005455SLois Curfman McInnes Notes: 1568c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1569c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1570c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1571c9005455SLois Curfman McInnes at each step. 1572c9005455SLois Curfman McInnes 1573c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1574c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1575c9005455SLois Curfman McInnes during the section of code that is being timed. 1576c9005455SLois Curfman McInnes 1577c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1578c9005455SLois Curfman McInnes @*/ 1579c9005455SLois Curfman McInnes int SNESSetConvergenceHistory(SNES snes, double *a, int na) 1580c9005455SLois Curfman McInnes { 15813a40ed3dSBarry Smith PetscFunctionBegin; 1582c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1583c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1584c9005455SLois Curfman McInnes snes->conv_hist = a; 1585c9005455SLois Curfman McInnes snes->conv_hist_size = na; 15863a40ed3dSBarry Smith PetscFunctionReturn(0); 1587c9005455SLois Curfman McInnes } 1588c9005455SLois Curfman McInnes 1589c9005455SLois Curfman McInnes #undef __FUNC__ 15905615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 15919b94acceSBarry Smith /* 15929b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 15939b94acceSBarry Smith positive parameter delta. 15949b94acceSBarry Smith 15959b94acceSBarry Smith Input Parameters: 1596c7afd0dbSLois Curfman McInnes + snes - the SNES context 15979b94acceSBarry Smith . y - approximate solution of linear system 15989b94acceSBarry Smith . fnorm - 2-norm of current function 1599c7afd0dbSLois Curfman McInnes - delta - trust region size 16009b94acceSBarry Smith 16019b94acceSBarry Smith Output Parameters: 1602c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 16039b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 16049b94acceSBarry Smith region, and exceeds zero otherwise. 1605c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 16069b94acceSBarry Smith 16079b94acceSBarry Smith Note: 160840191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 16099b94acceSBarry Smith is set to be the maximum allowable step size. 16109b94acceSBarry Smith 16119b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 16129b94acceSBarry Smith */ 16139b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 16149b94acceSBarry Smith double *gpnorm,double *ynorm) 16159b94acceSBarry Smith { 16169b94acceSBarry Smith double norm; 16179b94acceSBarry Smith Scalar cnorm; 16183a40ed3dSBarry Smith int ierr; 16193a40ed3dSBarry Smith 16203a40ed3dSBarry Smith PetscFunctionBegin; 16213a40ed3dSBarry Smith ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr); 16229b94acceSBarry Smith if (norm > *delta) { 16239b94acceSBarry Smith norm = *delta/norm; 16249b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 16259b94acceSBarry Smith cnorm = norm; 16269b94acceSBarry Smith VecScale( &cnorm, y ); 16279b94acceSBarry Smith *ynorm = *delta; 16289b94acceSBarry Smith } else { 16299b94acceSBarry Smith *gpnorm = 0.0; 16309b94acceSBarry Smith *ynorm = norm; 16319b94acceSBarry Smith } 16323a40ed3dSBarry Smith PetscFunctionReturn(0); 16339b94acceSBarry Smith } 16349b94acceSBarry Smith 16355615d1e5SSatish Balay #undef __FUNC__ 16365615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 16379b94acceSBarry Smith /*@ 16389b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 163963a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 16409b94acceSBarry Smith 1641c7afd0dbSLois Curfman McInnes Collective on SNES 1642c7afd0dbSLois Curfman McInnes 1643b2002411SLois Curfman McInnes Input Parameters: 1644c7afd0dbSLois Curfman McInnes + snes - the SNES context 1645c7afd0dbSLois Curfman McInnes - x - the solution vector 16469b94acceSBarry Smith 16479b94acceSBarry Smith Output Parameter: 1648b2002411SLois Curfman McInnes . its - number of iterations until termination 16499b94acceSBarry Smith 1650b2002411SLois Curfman McInnes Notes: 16518ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 16528ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 16538ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 16548ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 16558ddd3da0SLois Curfman McInnes 16569b94acceSBarry Smith .keywords: SNES, nonlinear, solve 16579b94acceSBarry Smith 165863a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 16599b94acceSBarry Smith @*/ 16608ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 16619b94acceSBarry Smith { 16623c7409f5SSatish Balay int ierr, flg; 1663052efed2SBarry Smith 16643a40ed3dSBarry Smith PetscFunctionBegin; 166577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 166674679c65SBarry Smith PetscValidIntPointer(its); 166782bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1668c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 16699b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 1670c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 16719b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 16729b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 16733c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 16746d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 16753a40ed3dSBarry Smith PetscFunctionReturn(0); 16769b94acceSBarry Smith } 16779b94acceSBarry Smith 16789b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 16799b94acceSBarry Smith 16805615d1e5SSatish Balay #undef __FUNC__ 16815615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 168282bf6240SBarry Smith /*@C 16834b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 16849b94acceSBarry Smith 1685fee21e36SBarry Smith Collective on SNES 1686fee21e36SBarry Smith 1687c7afd0dbSLois Curfman McInnes Input Parameters: 1688c7afd0dbSLois Curfman McInnes + snes - the SNES context 1689c7afd0dbSLois Curfman McInnes - method - a known method 1690c7afd0dbSLois Curfman McInnes 1691c7afd0dbSLois Curfman McInnes Options Database Key: 1692c7afd0dbSLois Curfman McInnes . -snes_type <method> - Sets the method; use -help for a list 1693c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1694ae12b187SLois Curfman McInnes 16959b94acceSBarry Smith Notes: 16969b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 1697c7afd0dbSLois Curfman McInnes . SNES_EQ_LS - Newton's method with line search 1698c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 1699c7afd0dbSLois Curfman McInnes . SNES_EQ_TR - Newton's method with trust region 1700c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 1701c7afd0dbSLois Curfman McInnes . SNES_UM_TR - Newton's method with trust region 1702c7afd0dbSLois Curfman McInnes (unconstrained minimization) 1703c7afd0dbSLois Curfman McInnes . SNES_UM_LS - Newton's method with line search 1704c7afd0dbSLois Curfman McInnes (unconstrained minimization) 17059b94acceSBarry Smith 1706ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1707ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1708ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1709ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1710ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1711ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1712ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1713ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1714ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 1715ae12b187SLois Curfman McInnes appropriate method. In other words, this routine is for the advanced user. 1716a703fe33SLois Curfman McInnes 1717f536c699SSatish Balay .keywords: SNES, set, method 17189b94acceSBarry Smith @*/ 17194b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 17209b94acceSBarry Smith { 172184cb2905SBarry Smith int ierr; 17229b94acceSBarry Smith int (*r)(SNES); 17233a40ed3dSBarry Smith 17243a40ed3dSBarry Smith PetscFunctionBegin; 172577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 172682bf6240SBarry Smith 172782bf6240SBarry Smith if (!PetscStrcmp(snes->type_name,method)) PetscFunctionReturn(0); 172882bf6240SBarry Smith 172982bf6240SBarry Smith if (snes->setupcalled) { 1730e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 173182bf6240SBarry Smith snes->data = 0; 173282bf6240SBarry Smith } 17337d1a2b2bSBarry Smith 17349b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 173582bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);} 173682bf6240SBarry Smith 1737ecf371e4SBarry Smith ierr = DLRegisterFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr); 173882bf6240SBarry Smith 17390452661fSBarry Smith if (snes->data) PetscFree(snes->data); 174082bf6240SBarry Smith snes->data = 0; 17413a40ed3dSBarry Smith ierr = (*r)(snes); CHKERRQ(ierr); 174282bf6240SBarry Smith 174382bf6240SBarry Smith if (snes->type_name) PetscFree(snes->type_name); 174482bf6240SBarry Smith snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name); 174582bf6240SBarry Smith PetscStrcpy(snes->type_name,method); 174682bf6240SBarry Smith snes->set_method_called = 1; 174782bf6240SBarry Smith 17483a40ed3dSBarry Smith PetscFunctionReturn(0); 17499b94acceSBarry Smith } 17509b94acceSBarry Smith 1751a847f771SSatish Balay 17529b94acceSBarry Smith /* --------------------------------------------------------------------- */ 17535615d1e5SSatish Balay #undef __FUNC__ 1754d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy" 17559b94acceSBarry Smith /*@C 17569b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 17579b94acceSBarry Smith registered by SNESRegister(). 17589b94acceSBarry Smith 1759fee21e36SBarry Smith Not Collective 1760fee21e36SBarry Smith 17619b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 17629b94acceSBarry Smith 17639b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 17649b94acceSBarry Smith @*/ 1765cf256101SBarry Smith int SNESRegisterDestroy(void) 17669b94acceSBarry Smith { 176782bf6240SBarry Smith int ierr; 176882bf6240SBarry Smith 17693a40ed3dSBarry Smith PetscFunctionBegin; 177082bf6240SBarry Smith if (SNESList) { 177182bf6240SBarry Smith ierr = DLRegisterDestroy( SNESList );CHKERRQ(ierr); 177282bf6240SBarry Smith SNESList = 0; 17739b94acceSBarry Smith } 177484cb2905SBarry Smith SNESRegisterAllCalled = 0; 17753a40ed3dSBarry Smith PetscFunctionReturn(0); 17769b94acceSBarry Smith } 17779b94acceSBarry Smith 17785615d1e5SSatish Balay #undef __FUNC__ 1779d4bb536fSBarry Smith #define __FUNC__ "SNESGetType" 17809b94acceSBarry Smith /*@C 17819a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 17829b94acceSBarry Smith 1783c7afd0dbSLois Curfman McInnes Not Collective 1784c7afd0dbSLois Curfman McInnes 17859b94acceSBarry Smith Input Parameter: 17864b0e389bSBarry Smith . snes - nonlinear solver context 17879b94acceSBarry Smith 17889b94acceSBarry Smith Output Parameter: 178982bf6240SBarry Smith . method - SNES method (a charactor string) 17909b94acceSBarry Smith 17919b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 17929b94acceSBarry Smith @*/ 179382bf6240SBarry Smith int SNESGetType(SNES snes, SNESType *method) 17949b94acceSBarry Smith { 17953a40ed3dSBarry Smith PetscFunctionBegin; 179682bf6240SBarry Smith *method = snes->type_name; 17973a40ed3dSBarry Smith PetscFunctionReturn(0); 17989b94acceSBarry Smith } 17999b94acceSBarry Smith 18005615d1e5SSatish Balay #undef __FUNC__ 1801d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution" 18029b94acceSBarry Smith /*@C 18039b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 18049b94acceSBarry Smith stored. 18059b94acceSBarry Smith 1806c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1807c7afd0dbSLois Curfman McInnes 18089b94acceSBarry Smith Input Parameter: 18099b94acceSBarry Smith . snes - the SNES context 18109b94acceSBarry Smith 18119b94acceSBarry Smith Output Parameter: 18129b94acceSBarry Smith . x - the solution 18139b94acceSBarry Smith 18149b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 18159b94acceSBarry Smith 1816a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 18179b94acceSBarry Smith @*/ 18189b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 18199b94acceSBarry Smith { 18203a40ed3dSBarry Smith PetscFunctionBegin; 182177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18229b94acceSBarry Smith *x = snes->vec_sol_always; 18233a40ed3dSBarry Smith PetscFunctionReturn(0); 18249b94acceSBarry Smith } 18259b94acceSBarry Smith 18265615d1e5SSatish Balay #undef __FUNC__ 1827d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" 18289b94acceSBarry Smith /*@C 18299b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 18309b94acceSBarry Smith stored. 18319b94acceSBarry Smith 1832c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1833c7afd0dbSLois Curfman McInnes 18349b94acceSBarry Smith Input Parameter: 18359b94acceSBarry Smith . snes - the SNES context 18369b94acceSBarry Smith 18379b94acceSBarry Smith Output Parameter: 18389b94acceSBarry Smith . x - the solution update 18399b94acceSBarry Smith 18409b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 18419b94acceSBarry Smith 18429b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 18439b94acceSBarry Smith @*/ 18449b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 18459b94acceSBarry Smith { 18463a40ed3dSBarry Smith PetscFunctionBegin; 184777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18489b94acceSBarry Smith *x = snes->vec_sol_update_always; 18493a40ed3dSBarry Smith PetscFunctionReturn(0); 18509b94acceSBarry Smith } 18519b94acceSBarry Smith 18525615d1e5SSatish Balay #undef __FUNC__ 1853d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction" 18549b94acceSBarry Smith /*@C 18553638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 18569b94acceSBarry Smith 1857c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1858c7afd0dbSLois Curfman McInnes 18599b94acceSBarry Smith Input Parameter: 18609b94acceSBarry Smith . snes - the SNES context 18619b94acceSBarry Smith 18629b94acceSBarry Smith Output Parameter: 18633638b69dSLois Curfman McInnes . r - the function 18649b94acceSBarry Smith 18659b94acceSBarry Smith Notes: 18669b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 18679b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 18689b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 18699b94acceSBarry Smith 1870a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 18719b94acceSBarry Smith 187231615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 187331615d04SLois Curfman McInnes SNESGetGradient() 18749b94acceSBarry Smith @*/ 18759b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 18769b94acceSBarry Smith { 18773a40ed3dSBarry Smith PetscFunctionBegin; 187877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1879a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1880a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1881a8c6a408SBarry Smith } 18829b94acceSBarry Smith *r = snes->vec_func_always; 18833a40ed3dSBarry Smith PetscFunctionReturn(0); 18849b94acceSBarry Smith } 18859b94acceSBarry Smith 18865615d1e5SSatish Balay #undef __FUNC__ 1887d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient" 18889b94acceSBarry Smith /*@C 18893638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 18909b94acceSBarry Smith 1891c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 1892c7afd0dbSLois Curfman McInnes 18939b94acceSBarry Smith Input Parameter: 18949b94acceSBarry Smith . snes - the SNES context 18959b94acceSBarry Smith 18969b94acceSBarry Smith Output Parameter: 18979b94acceSBarry Smith . r - the gradient 18989b94acceSBarry Smith 18999b94acceSBarry Smith Notes: 19009b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 19019b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 19029b94acceSBarry Smith SNESGetFunction(). 19039b94acceSBarry Smith 19049b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 19059b94acceSBarry Smith 190631615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 19079b94acceSBarry Smith @*/ 19089b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 19099b94acceSBarry Smith { 19103a40ed3dSBarry Smith PetscFunctionBegin; 191177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 19123a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1913a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 19143a40ed3dSBarry Smith } 19159b94acceSBarry Smith *r = snes->vec_func_always; 19163a40ed3dSBarry Smith PetscFunctionReturn(0); 19179b94acceSBarry Smith } 19189b94acceSBarry Smith 19195615d1e5SSatish Balay #undef __FUNC__ 1920d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" 19219b94acceSBarry Smith /*@ 19229b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 19239b94acceSBarry Smith unconstrained minimization problems. 19249b94acceSBarry Smith 1925c7afd0dbSLois Curfman McInnes Not Collective 1926c7afd0dbSLois Curfman McInnes 19279b94acceSBarry Smith Input Parameter: 19289b94acceSBarry Smith . snes - the SNES context 19299b94acceSBarry Smith 19309b94acceSBarry Smith Output Parameter: 19319b94acceSBarry Smith . r - the function 19329b94acceSBarry Smith 19339b94acceSBarry Smith Notes: 19349b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 19359b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 19369b94acceSBarry Smith SNESGetFunction(). 19379b94acceSBarry Smith 19389b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 19399b94acceSBarry Smith 194031615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 19419b94acceSBarry Smith @*/ 19429b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 19439b94acceSBarry Smith { 19443a40ed3dSBarry Smith PetscFunctionBegin; 194577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 194674679c65SBarry Smith PetscValidScalarPointer(r); 19473a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1948a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 19493a40ed3dSBarry Smith } 19509b94acceSBarry Smith *r = snes->fc; 19513a40ed3dSBarry Smith PetscFunctionReturn(0); 19529b94acceSBarry Smith } 19539b94acceSBarry Smith 19545615d1e5SSatish Balay #undef __FUNC__ 1955d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" 19563c7409f5SSatish Balay /*@C 19573c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 1958d850072dSLois Curfman McInnes SNES options in the database. 19593c7409f5SSatish Balay 1960fee21e36SBarry Smith Collective on SNES 1961fee21e36SBarry Smith 1962c7afd0dbSLois Curfman McInnes Input Parameter: 1963c7afd0dbSLois Curfman McInnes + snes - the SNES context 1964c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1965c7afd0dbSLois Curfman McInnes 1966d850072dSLois Curfman McInnes Notes: 1967a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1968c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 1969d850072dSLois Curfman McInnes 19703c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 1971a86d99e1SLois Curfman McInnes 1972a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 19733c7409f5SSatish Balay @*/ 19743c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 19753c7409f5SSatish Balay { 19763c7409f5SSatish Balay int ierr; 19773c7409f5SSatish Balay 19783a40ed3dSBarry Smith PetscFunctionBegin; 197977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1980639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 19813c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 19823a40ed3dSBarry Smith PetscFunctionReturn(0); 19833c7409f5SSatish Balay } 19843c7409f5SSatish Balay 19855615d1e5SSatish Balay #undef __FUNC__ 1986d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" 19873c7409f5SSatish Balay /*@C 1988f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 1989d850072dSLois Curfman McInnes SNES options in the database. 19903c7409f5SSatish Balay 1991fee21e36SBarry Smith Collective on SNES 1992fee21e36SBarry Smith 1993c7afd0dbSLois Curfman McInnes Input Parameters: 1994c7afd0dbSLois Curfman McInnes + snes - the SNES context 1995c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 1996c7afd0dbSLois Curfman McInnes 1997d850072dSLois Curfman McInnes Notes: 1998a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 1999c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2000d850072dSLois Curfman McInnes 20013c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2002a86d99e1SLois Curfman McInnes 2003a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 20043c7409f5SSatish Balay @*/ 20053c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 20063c7409f5SSatish Balay { 20073c7409f5SSatish Balay int ierr; 20083c7409f5SSatish Balay 20093a40ed3dSBarry Smith PetscFunctionBegin; 201077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2011639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 20123c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 20133a40ed3dSBarry Smith PetscFunctionReturn(0); 20143c7409f5SSatish Balay } 20153c7409f5SSatish Balay 20165615d1e5SSatish Balay #undef __FUNC__ 2017d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" 2018c83590e2SSatish Balay /*@ 20193c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 20203c7409f5SSatish Balay SNES options in the database. 20213c7409f5SSatish Balay 2022c7afd0dbSLois Curfman McInnes Not Collective 2023c7afd0dbSLois Curfman McInnes 20243c7409f5SSatish Balay Input Parameter: 20253c7409f5SSatish Balay . snes - the SNES context 20263c7409f5SSatish Balay 20273c7409f5SSatish Balay Output Parameter: 20283c7409f5SSatish Balay . prefix - pointer to the prefix string used 20293c7409f5SSatish Balay 20303c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2031a86d99e1SLois Curfman McInnes 2032a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 20333c7409f5SSatish Balay @*/ 20343c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 20353c7409f5SSatish Balay { 20363c7409f5SSatish Balay int ierr; 20373c7409f5SSatish Balay 20383a40ed3dSBarry Smith PetscFunctionBegin; 203977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2040639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 20413a40ed3dSBarry Smith PetscFunctionReturn(0); 20423c7409f5SSatish Balay } 20433c7409f5SSatish Balay 2044ca161407SBarry Smith #undef __FUNC__ 2045ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp" 2046ca161407SBarry Smith /*@ 2047ca161407SBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2048ca161407SBarry Smith 2049c7afd0dbSLois Curfman McInnes Collective on SNES 2050c7afd0dbSLois Curfman McInnes 2051ca161407SBarry Smith Input Parameter: 2052ca161407SBarry Smith . snes - the SNES context 2053ca161407SBarry Smith 2054ca161407SBarry Smith Options Database Keys: 2055c7afd0dbSLois Curfman McInnes + -help - Prints SNES options 2056c7afd0dbSLois Curfman McInnes - -h - Prints SNES options 2057ca161407SBarry Smith 2058ca161407SBarry Smith .keywords: SNES, nonlinear, help 2059ca161407SBarry Smith 2060ca161407SBarry Smith .seealso: SNESSetFromOptions() 2061ca161407SBarry Smith @*/ 2062ca161407SBarry Smith int SNESPrintHelp(SNES snes) 2063ca161407SBarry Smith { 2064ca161407SBarry Smith char p[64]; 2065ca161407SBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2066ca161407SBarry Smith int ierr; 2067ca161407SBarry Smith 2068ca161407SBarry Smith PetscFunctionBegin; 2069ca161407SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2070ca161407SBarry Smith 2071ca161407SBarry Smith PetscStrcpy(p,"-"); 2072ca161407SBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2073ca161407SBarry Smith 2074ca161407SBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2075ca161407SBarry Smith 207676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n"); 207782bf6240SBarry Smith ierr = DLRegisterPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr); 207876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 207976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 208076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 208176be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 208276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 208376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 208476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 208576be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 208676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 208776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 2088ca161407SBarry Smith residual norm at each iteration.\n",p); 208976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 2090ca161407SBarry Smith residual norm for small residual norms. This is useful to conceal\n\ 2091ca161407SBarry Smith meaningless digits that may be different on different machines.\n",p); 209276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 2093ca161407SBarry Smith if (snes->type == SNES_NONLINEAR_EQUATIONS) { 209476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2095ca161407SBarry Smith " Options for solving systems of nonlinear equations only:\n"); 209676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 209776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 209876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 209976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p); 210076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 210176be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2102ca161407SBarry Smith " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 210376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2104ca161407SBarry Smith " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 210576be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2106ca161407SBarry Smith " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 210776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2108ca161407SBarry Smith " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 210976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2110ca161407SBarry Smith " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 211176be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2112ca161407SBarry Smith " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 211376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2114ca161407SBarry Smith " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 2115ca161407SBarry Smith } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 211676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2117ca161407SBarry Smith " Options for solving unconstrained minimization problems only:\n"); 211876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 211976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 212076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 2121ca161407SBarry Smith } 212276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 212376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm,"a particular method\n"); 2124ca161407SBarry Smith if (snes->printhelp) { 2125ca161407SBarry Smith ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr); 2126ca161407SBarry Smith } 2127ca161407SBarry Smith PetscFunctionReturn(0); 2128ca161407SBarry Smith } 21293c7409f5SSatish Balay 2130acb85ae6SSatish Balay /*MC 2131b2002411SLois Curfman McInnes SNESRegister - Adds a method to the nonlinear solver package. 21329b94acceSBarry Smith 2133b2002411SLois Curfman McInnes Synopsis: 2134b2002411SLois Curfman McInnes SNESRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES)) 21359b94acceSBarry Smith 2136b2002411SLois Curfman McInnes Input Parameters: 2137c7afd0dbSLois Curfman McInnes + name_solver - name of a new user-defined solver 2138b2002411SLois Curfman McInnes . path - path (either absolute or relative) the library containing this solver 2139b2002411SLois Curfman McInnes . name_create - name of routine to create method context 2140c7afd0dbSLois Curfman McInnes - routine_create - routine to create method context 21419b94acceSBarry Smith 2142b2002411SLois Curfman McInnes Notes: 2143b2002411SLois Curfman McInnes SNESRegister() may be called multiple times to add several user-defined solvers. 21443a40ed3dSBarry Smith 2145b2002411SLois Curfman McInnes If dynamic libraries are used, then the fourth input argument (routine_create) 2146b2002411SLois Curfman McInnes is ignored. 2147b2002411SLois Curfman McInnes 2148b2002411SLois Curfman McInnes Sample usage: 214969225d78SLois Curfman McInnes .vb 2150b2002411SLois Curfman McInnes SNESRegister("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a, 2151b2002411SLois Curfman McInnes "MySolverCreate",MySolverCreate); 215269225d78SLois Curfman McInnes .ve 2153b2002411SLois Curfman McInnes 2154b2002411SLois Curfman McInnes Then, your solver can be chosen with the procedural interface via 2155b2002411SLois Curfman McInnes $ SNESSetType(snes,"my_solver") 2156b2002411SLois Curfman McInnes or at runtime via the option 2157b2002411SLois Curfman McInnes $ -snes_type my_solver 2158b2002411SLois Curfman McInnes 2159b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register 2160b2002411SLois Curfman McInnes 2161b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy() 2162b2002411SLois Curfman McInnes M*/ 2163b2002411SLois Curfman McInnes 2164b2002411SLois Curfman McInnes #undef __FUNC__ 2165b2002411SLois Curfman McInnes #define __FUNC__ "SNESRegister_Private" 2166b2002411SLois Curfman McInnes int SNESRegister_Private(char *sname,char *path,char *name,int (*function)(SNES)) 2167b2002411SLois Curfman McInnes { 2168b2002411SLois Curfman McInnes char fullname[256]; 2169b2002411SLois Curfman McInnes int ierr; 2170b2002411SLois Curfman McInnes 2171b2002411SLois Curfman McInnes PetscFunctionBegin; 2172b2002411SLois Curfman McInnes PetscStrcpy(fullname,path); PetscStrcat(fullname,":");PetscStrcat(fullname,name); 2173b2002411SLois Curfman McInnes ierr = DLRegister_Private(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr); 2174b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2175b2002411SLois Curfman McInnes } 2176