1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*9ab63eb5SSatish Balay static char vcid[] = "$Id: snes.c,v 1.181 1999/04/05 18:18:42 bsmith Exp balay $"; 39b94acceSBarry Smith #endif 49b94acceSBarry Smith 570f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 69b94acceSBarry Smith 784cb2905SBarry Smith int SNESRegisterAllCalled = 0; 8488ecbafSBarry Smith FList SNESList = 0; 982bf6240SBarry Smith 105615d1e5SSatish Balay #undef __FUNC__ 11d4bb536fSBarry Smith #define __FUNC__ "SNESView" 129b94acceSBarry Smith /*@ 139b94acceSBarry Smith SNESView - Prints the SNES data structure. 149b94acceSBarry Smith 15fee21e36SBarry Smith Collective on SNES, unless Viewer is VIEWER_STDOUT_SELF 16fee21e36SBarry Smith 17c7afd0dbSLois Curfman McInnes Input Parameters: 18c7afd0dbSLois Curfman McInnes + SNES - the SNES context 19c7afd0dbSLois Curfman McInnes - viewer - visualization context 20c7afd0dbSLois Curfman McInnes 219b94acceSBarry Smith Options Database Key: 22c8a8ba5cSLois Curfman McInnes . -snes_view - Calls SNESView() at end of SNESSolve() 239b94acceSBarry Smith 249b94acceSBarry Smith Notes: 259b94acceSBarry Smith The available visualization contexts include 26c7afd0dbSLois Curfman McInnes + VIEWER_STDOUT_SELF - standard output (default) 27c7afd0dbSLois Curfman McInnes - VIEWER_STDOUT_WORLD - synchronized standard 28c8a8ba5cSLois Curfman McInnes output where only the first processor opens 29c8a8ba5cSLois Curfman McInnes the file. All other processors send their 30c8a8ba5cSLois Curfman McInnes data to the first processor to print. 319b94acceSBarry Smith 323e081fefSLois Curfman McInnes The user can open an alternative visualization context with 3377ed5343SBarry Smith ViewerASCIIOpen() - output to a specified file. 349b94acceSBarry Smith 3536851e7fSLois Curfman McInnes Level: beginner 3636851e7fSLois Curfman McInnes 379b94acceSBarry Smith .keywords: SNES, view 389b94acceSBarry Smith 3977ed5343SBarry Smith .seealso: ViewerASCIIOpen() 409b94acceSBarry Smith @*/ 419b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 429b94acceSBarry Smith { 439b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 449b94acceSBarry Smith int ierr; 459b94acceSBarry Smith SLES sles; 469b94acceSBarry Smith char *method; 4719bcc07fSBarry Smith ViewerType vtype; 489b94acceSBarry Smith 493a40ed3dSBarry Smith PetscFunctionBegin; 5074679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5174679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5274679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5374679c65SBarry Smith 5419bcc07fSBarry Smith ierr = ViewerGetType(viewer,&vtype); CHKERRQ(ierr); 553f1db9ecSBarry Smith if (PetscTypeCompare(vtype,ASCII_VIEWER)) { 560ef38995SBarry Smith ViewerASCIIPrintf(viewer,"SNES Object:\n"); 5782bf6240SBarry Smith SNESGetType(snes,&method); 580ef38995SBarry Smith ViewerASCIIPrintf(viewer," method: %s\n",method); 590ef38995SBarry Smith if (snes->view) { 600ef38995SBarry Smith ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 610ef38995SBarry Smith ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 620ef38995SBarry Smith ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 630ef38995SBarry Smith } 640ef38995SBarry Smith ViewerASCIIPrintf(viewer," maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs); 650ef38995SBarry Smith ViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 669b94acceSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol); 670ef38995SBarry Smith ViewerASCIIPrintf(viewer," total number of linear solver iterations=%d\n",snes->linear_its); 680ef38995SBarry Smith ViewerASCIIPrintf(viewer," total number of function evaluations=%d\n",snes->nfuncs); 690ef38995SBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 700ef38995SBarry Smith ViewerASCIIPrintf(viewer," min function tolerance=%g\n",snes->fmin); 710ef38995SBarry Smith } 729b94acceSBarry Smith if (snes->ksp_ewconv) { 739b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 749b94acceSBarry Smith if (kctx) { 750ef38995SBarry Smith ViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version); 760ef38995SBarry Smith ViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold); 770ef38995SBarry Smith ViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2); 789b94acceSBarry Smith } 799b94acceSBarry Smith } 803f1db9ecSBarry Smith } else if (PetscTypeCompare(vtype,STRING_VIEWER)) { 810ef38995SBarry Smith ierr = SNESGetType(snes,&method);CHKERRQ(ierr); 820ef38995SBarry Smith ierr = ViewerStringSPrintf(viewer," %-3.3s",method);CHKERRQ(ierr); 8319bcc07fSBarry Smith } 8477ed5343SBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 850ef38995SBarry Smith ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 869b94acceSBarry Smith ierr = SLESView(sles,viewer); CHKERRQ(ierr); 870ef38995SBarry Smith ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 883a40ed3dSBarry Smith PetscFunctionReturn(0); 899b94acceSBarry Smith } 909b94acceSBarry Smith 91639f9d9dSBarry Smith /* 92639f9d9dSBarry Smith We retain a list of functions that also take SNES command 93639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 94639f9d9dSBarry Smith */ 95639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 96639f9d9dSBarry Smith static int numberofsetfromoptions; 97639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 98639f9d9dSBarry Smith 995615d1e5SSatish Balay #undef __FUNC__ 100d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker" 101639f9d9dSBarry Smith /*@ 102639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 103639f9d9dSBarry Smith 104c7afd0dbSLois Curfman McInnes Not Collective 105c7afd0dbSLois Curfman McInnes 106639f9d9dSBarry Smith Input Parameter: 107639f9d9dSBarry Smith . snescheck - function that checks for options 108639f9d9dSBarry Smith 10936851e7fSLois Curfman McInnes Level: developer 11036851e7fSLois Curfman McInnes 111639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 112639f9d9dSBarry Smith @*/ 113639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 114639f9d9dSBarry Smith { 1153a40ed3dSBarry Smith PetscFunctionBegin; 116639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 117a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed"); 118639f9d9dSBarry Smith } 119639f9d9dSBarry Smith 120639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 1213a40ed3dSBarry Smith PetscFunctionReturn(0); 122639f9d9dSBarry Smith } 123639f9d9dSBarry Smith 1245615d1e5SSatish Balay #undef __FUNC__ 12515091d37SBarry Smith #define __FUNC__ "SNESSetTypeFromOptions" 12615091d37SBarry Smith /*@ 12715091d37SBarry Smith SNESSetTypeFromOptions - Sets the SNES solver type from the options database, 12815091d37SBarry Smith or sets a default if none is give. 12915091d37SBarry Smith 13015091d37SBarry Smith Collective on SNES 13115091d37SBarry Smith 13215091d37SBarry Smith Input Parameter: 13315091d37SBarry Smith . snes - the SNES context 13415091d37SBarry Smith 13515091d37SBarry Smith Options Database Keys: 13615091d37SBarry Smith . -snes_type <type> - SNES_EQ_LS, SNES_EQ_TR, SNES_UM_TR, SNES_UM_LS etc 13715091d37SBarry Smith 13815091d37SBarry Smith Level: beginner 13915091d37SBarry Smith 14015091d37SBarry Smith .keywords: SNES, nonlinear, set, options, database 14115091d37SBarry Smith 14215091d37SBarry Smith .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetFromOptions() 14315091d37SBarry Smith @*/ 14415091d37SBarry Smith int SNESSetTypeFromOptions(SNES snes) 14515091d37SBarry Smith { 14615091d37SBarry Smith char method[256]; 14715091d37SBarry Smith int ierr, flg; 14815091d37SBarry Smith 14915091d37SBarry Smith PetscFunctionBegin; 15015091d37SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 15115091d37SBarry Smith if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp()"); 15215091d37SBarry Smith ierr = OptionsGetString(snes->prefix,"-snes_type",method,256,&flg); 15315091d37SBarry Smith if (flg) { 15415091d37SBarry Smith ierr = SNESSetType(snes,(SNESType) method); CHKERRQ(ierr); 15515091d37SBarry Smith } 15615091d37SBarry Smith /* 15715091d37SBarry Smith If SNES type has not yet been set, set it now 15815091d37SBarry Smith */ 15915091d37SBarry Smith if (!snes->type_name) { 16015091d37SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 16115091d37SBarry Smith ierr = SNESSetType(snes,SNES_EQ_LS); CHKERRQ(ierr); 16215091d37SBarry Smith } else { 16315091d37SBarry Smith ierr = SNESSetType(snes,SNES_UM_TR); CHKERRQ(ierr); 16415091d37SBarry Smith } 16515091d37SBarry Smith } 16615091d37SBarry Smith PetscFunctionReturn(0); 16715091d37SBarry Smith } 16815091d37SBarry Smith 16915091d37SBarry Smith #undef __FUNC__ 1705615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1719b94acceSBarry Smith /*@ 172682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1739b94acceSBarry Smith 174c7afd0dbSLois Curfman McInnes Collective on SNES 175c7afd0dbSLois Curfman McInnes 1769b94acceSBarry Smith Input Parameter: 1779b94acceSBarry Smith . snes - the SNES context 1789b94acceSBarry Smith 17936851e7fSLois Curfman McInnes Options Database Keys: 180b39c3a46SLois Curfman McInnes + -snes_type <type> - SNES_EQ_LS, SNES_EQ_TR, SNES_UM_TR, SNES_UM_LS etc 18182738288SBarry Smith . -snes_stol - convergence tolerance in terms of the norm 18282738288SBarry Smith of the change in the solution between steps 183b39c3a46SLois Curfman McInnes . -snes_atol <atol> - absolute tolerance of residual norm 184b39c3a46SLois Curfman McInnes . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 185b39c3a46SLois Curfman McInnes . -snes_max_it <max_it> - maximum number of iterations 186b39c3a46SLois Curfman McInnes . -snes_max_funcs <max_funcs> - maximum number of function evaluations 187b39c3a46SLois Curfman McInnes . -snes_trtol <trtol> - trust region tolerance 18882738288SBarry Smith . -snes_no_convergence_test - skip convergence test in nonlinear or minimization 18982738288SBarry Smith solver; hence iterations will continue until max_it 1901fbbfb26SLois Curfman McInnes or some other criterion is reached. Saves expense 19182738288SBarry Smith of convergence test 19282738288SBarry Smith . -snes_monitor - prints residual norm at each iteration 19382738288SBarry Smith . -snes_xmonitor - plots residual norm at each iteration 194e24b481bSBarry Smith . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 19536851e7fSLois Curfman McInnes - -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 19682738288SBarry Smith 19782738288SBarry Smith Options Database for Eisenstat-Walker method: 19882738288SBarry Smith + -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence 19982738288SBarry Smith . -snes_ksp_eq_version ver - version of Eisenstat-Walker method 20036851e7fSLois Curfman McInnes . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 20136851e7fSLois Curfman McInnes . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 20236851e7fSLois Curfman McInnes . -snes_ksp_ew_gamma <gamma> - Sets gamma 20336851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha <alpha> - Sets alpha 20436851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 20536851e7fSLois Curfman McInnes - -snes_ksp_ew_threshold <threshold> - Sets threshold 20682738288SBarry Smith 20711ca99fdSLois Curfman McInnes Notes: 20811ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 20911ca99fdSLois Curfman McInnes the users manual. 21083e2fdc7SBarry Smith 21136851e7fSLois Curfman McInnes Level: beginner 21236851e7fSLois Curfman McInnes 2139b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 2149b94acceSBarry Smith 21515091d37SBarry Smith .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetTypeFromOptions() 2169b94acceSBarry Smith @*/ 2179b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 2189b94acceSBarry Smith { 2199b94acceSBarry Smith double tmp; 2209b94acceSBarry Smith SLES sles; 22181f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 2223c7409f5SSatish Balay int version = PETSC_DEFAULT; 2239b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 2249b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 2259b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 2269b94acceSBarry Smith double alpha = PETSC_DEFAULT; 2279b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 2289b94acceSBarry Smith double threshold = PETSC_DEFAULT; 2299b94acceSBarry Smith 2303a40ed3dSBarry Smith PetscFunctionBegin; 23177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 23215091d37SBarry Smith ierr = SNESSetTypeFromOptions(snes);CHKERRQ(ierr); 233ca161407SBarry Smith 23415091d37SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 2353c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg); CHKERRQ(ierr); 236d64ed03dSBarry Smith if (flg) { 237d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 238d64ed03dSBarry Smith } 2393c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg); CHKERRQ(ierr); 240d64ed03dSBarry Smith if (flg) { 241d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 242d64ed03dSBarry Smith } 2433c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg); CHKERRQ(ierr); 244d64ed03dSBarry Smith if (flg) { 245d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 246d64ed03dSBarry Smith } 2473c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg); CHKERRQ(ierr); 2483c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 249d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg); CHKERRQ(ierr); 250d64ed03dSBarry Smith if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp); CHKERRQ(ierr); } 251d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg); CHKERRQ(ierr); 252d64ed03dSBarry Smith if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp); CHKERRQ(ierr);} 2533c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg); CHKERRQ(ierr); 2543c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 255b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg); CHKERRQ(ierr); 256b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg); CHKERRQ(ierr); 257b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 258b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg); CHKERRQ(ierr); 259b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg); CHKERRQ(ierr); 260b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg); CHKERRQ(ierr); 261b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 26293c39befSBarry Smith 26393c39befSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg); CHKERRQ(ierr); 26493c39befSBarry Smith if (flg) {snes->converged = 0;} 26593c39befSBarry Smith 2669b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 2679b94acceSBarry Smith alpha2,threshold); CHKERRQ(ierr); 2685bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg); CHKERRQ(ierr); 2695cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 2703c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg); CHKERRQ(ierr); 271639f9d9dSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0);CHKERRQ(ierr);} 2723c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg); CHKERRQ(ierr); 2733f1db9ecSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0); CHKERRQ(ierr);} 2743f1db9ecSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor",&flg); CHKERRQ(ierr); 2753f1db9ecSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0); CHKERRQ(ierr);} 27681f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2773c7409f5SSatish Balay if (flg) { 27817699dbbSLois Curfman McInnes int rank = 0; 279d7e8b826SBarry Smith DrawLG lg; 28017699dbbSLois Curfman McInnes MPI_Initialized(&rank); 28117699dbbSLois Curfman McInnes if (rank) MPI_Comm_rank(snes->comm,&rank); 28217699dbbSLois Curfman McInnes if (!rank) { 28381f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg); CHKERRQ(ierr); 2849b94acceSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,(void *)lg); CHKERRQ(ierr); 285f8c826e1SBarry Smith snes->xmonitor = lg; 2869b94acceSBarry Smith PLogObjectParent(snes,lg); 2879b94acceSBarry Smith } 2889b94acceSBarry Smith } 289e24b481bSBarry Smith 2903c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg); CHKERRQ(ierr); 2913c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2929b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2939b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP); CHKERRQ(ierr); 294981c4779SBarry Smith PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 295981c4779SBarry Smith } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 29631615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 29731615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP); CHKERRQ(ierr); 298d64ed03dSBarry Smith PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 2999b94acceSBarry Smith } 300639f9d9dSBarry Smith 301639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 302639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes); CHKERRQ(ierr); 303639f9d9dSBarry Smith } 304639f9d9dSBarry Smith 3059b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 3069b94acceSBarry Smith ierr = SLESSetFromOptions(sles); CHKERRQ(ierr); 30793993e2dSLois Curfman McInnes 308e24b481bSBarry Smith /* set the special KSP monitor for matrix-free application */ 309e24b481bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_ksp_monitor",&flg); CHKERRQ(ierr); 310e24b481bSBarry Smith if (flg) { 311e24b481bSBarry Smith KSP ksp; 312e24b481bSBarry Smith ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 3135a655dc6SBarry Smith ierr = KSPSetMonitor(ksp,MatSNESMFKSPMonitor,PETSC_NULL);CHKERRQ(ierr); 314e24b481bSBarry Smith } 315e24b481bSBarry Smith 31682bf6240SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg); CHKERRQ(ierr); 31782bf6240SBarry Smith if (flg) { ierr = SNESPrintHelp(snes); CHKERRQ(ierr);} 31882bf6240SBarry Smith 3193a40ed3dSBarry Smith if (snes->setfromoptions) { 3203a40ed3dSBarry Smith ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 3213a40ed3dSBarry Smith } 3223a40ed3dSBarry Smith PetscFunctionReturn(0); 3239b94acceSBarry Smith } 3249b94acceSBarry Smith 325a847f771SSatish Balay 3265615d1e5SSatish Balay #undef __FUNC__ 327d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext" 3289b94acceSBarry Smith /*@ 3299b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 3309b94acceSBarry Smith the nonlinear solvers. 3319b94acceSBarry Smith 332fee21e36SBarry Smith Collective on SNES 333fee21e36SBarry Smith 334c7afd0dbSLois Curfman McInnes Input Parameters: 335c7afd0dbSLois Curfman McInnes + snes - the SNES context 336c7afd0dbSLois Curfman McInnes - usrP - optional user context 337c7afd0dbSLois Curfman McInnes 33836851e7fSLois Curfman McInnes Level: intermediate 33936851e7fSLois Curfman McInnes 3409b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3419b94acceSBarry Smith 3429b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3439b94acceSBarry Smith @*/ 3449b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3459b94acceSBarry Smith { 3463a40ed3dSBarry Smith PetscFunctionBegin; 34777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3489b94acceSBarry Smith snes->user = usrP; 3493a40ed3dSBarry Smith PetscFunctionReturn(0); 3509b94acceSBarry Smith } 35174679c65SBarry Smith 3525615d1e5SSatish Balay #undef __FUNC__ 353d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext" 3549b94acceSBarry Smith /*@C 3559b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3569b94acceSBarry Smith nonlinear solvers. 3579b94acceSBarry Smith 358c7afd0dbSLois Curfman McInnes Not Collective 359c7afd0dbSLois Curfman McInnes 3609b94acceSBarry Smith Input Parameter: 3619b94acceSBarry Smith . snes - SNES context 3629b94acceSBarry Smith 3639b94acceSBarry Smith Output Parameter: 3649b94acceSBarry Smith . usrP - user context 3659b94acceSBarry Smith 36636851e7fSLois Curfman McInnes Level: intermediate 36736851e7fSLois Curfman McInnes 3689b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3699b94acceSBarry Smith 3709b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3719b94acceSBarry Smith @*/ 3729b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3739b94acceSBarry Smith { 3743a40ed3dSBarry Smith PetscFunctionBegin; 37577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3769b94acceSBarry Smith *usrP = snes->user; 3773a40ed3dSBarry Smith PetscFunctionReturn(0); 3789b94acceSBarry Smith } 37974679c65SBarry Smith 3805615d1e5SSatish Balay #undef __FUNC__ 381d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber" 3829b94acceSBarry Smith /*@ 383c8228a4eSBarry Smith SNESGetIterationNumber - Gets the number of nonlinear iterations completed 384c8228a4eSBarry Smith at this time. 3859b94acceSBarry Smith 386c7afd0dbSLois Curfman McInnes Not Collective 387c7afd0dbSLois Curfman McInnes 3889b94acceSBarry Smith Input Parameter: 3899b94acceSBarry Smith . snes - SNES context 3909b94acceSBarry Smith 3919b94acceSBarry Smith Output Parameter: 3929b94acceSBarry Smith . iter - iteration number 3939b94acceSBarry Smith 394c8228a4eSBarry Smith Notes: 395c8228a4eSBarry Smith For example, during the computation of iteration 2 this would return 1. 396c8228a4eSBarry Smith 397c8228a4eSBarry Smith This is useful for using lagged Jacobians (where one does not recompute the 39808405cd6SLois Curfman McInnes Jacobian at each SNES iteration). For example, the code 39908405cd6SLois Curfman McInnes .vb 40008405cd6SLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&it); 40108405cd6SLois Curfman McInnes if (!(it % 2)) { 40208405cd6SLois Curfman McInnes [compute Jacobian here] 40308405cd6SLois Curfman McInnes } 40408405cd6SLois Curfman McInnes .ve 405c8228a4eSBarry Smith can be used in your ComputeJacobian() function to cause the Jacobian to be 40608405cd6SLois Curfman McInnes recomputed every second SNES iteration. 407c8228a4eSBarry Smith 40836851e7fSLois Curfman McInnes Level: intermediate 40936851e7fSLois Curfman McInnes 4109b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 4119b94acceSBarry Smith @*/ 4129b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 4139b94acceSBarry Smith { 4143a40ed3dSBarry Smith PetscFunctionBegin; 41577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 41674679c65SBarry Smith PetscValidIntPointer(iter); 4179b94acceSBarry Smith *iter = snes->iter; 4183a40ed3dSBarry Smith PetscFunctionReturn(0); 4199b94acceSBarry Smith } 42074679c65SBarry Smith 4215615d1e5SSatish Balay #undef __FUNC__ 4225615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 4239b94acceSBarry Smith /*@ 4249b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 4259b94acceSBarry Smith with SNESSSetFunction(). 4269b94acceSBarry Smith 427c7afd0dbSLois Curfman McInnes Collective on SNES 428c7afd0dbSLois Curfman McInnes 4299b94acceSBarry Smith Input Parameter: 4309b94acceSBarry Smith . snes - SNES context 4319b94acceSBarry Smith 4329b94acceSBarry Smith Output Parameter: 4339b94acceSBarry Smith . fnorm - 2-norm of function 4349b94acceSBarry Smith 4359b94acceSBarry Smith Note: 4369b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 437a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 438a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 4399b94acceSBarry Smith 44036851e7fSLois Curfman McInnes Level: intermediate 44136851e7fSLois Curfman McInnes 4429b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 443a86d99e1SLois Curfman McInnes 44408405cd6SLois Curfman McInnes .seealso: SNESGetFunction() 4459b94acceSBarry Smith @*/ 4469b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 4479b94acceSBarry Smith { 4483a40ed3dSBarry Smith PetscFunctionBegin; 44977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 45074679c65SBarry Smith PetscValidScalarPointer(fnorm); 45174679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 452d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 45374679c65SBarry Smith } 4549b94acceSBarry Smith *fnorm = snes->norm; 4553a40ed3dSBarry Smith PetscFunctionReturn(0); 4569b94acceSBarry Smith } 45774679c65SBarry Smith 4585615d1e5SSatish Balay #undef __FUNC__ 4595615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 4609b94acceSBarry Smith /*@ 4619b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 4629b94acceSBarry Smith with SNESSSetGradient(). 4639b94acceSBarry Smith 464c7afd0dbSLois Curfman McInnes Collective on SNES 465c7afd0dbSLois Curfman McInnes 4669b94acceSBarry Smith Input Parameter: 4679b94acceSBarry Smith . snes - SNES context 4689b94acceSBarry Smith 4699b94acceSBarry Smith Output Parameter: 4709b94acceSBarry Smith . fnorm - 2-norm of gradient 4719b94acceSBarry Smith 4729b94acceSBarry Smith Note: 4739b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 474a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 475a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4769b94acceSBarry Smith 47736851e7fSLois Curfman McInnes Level: intermediate 47836851e7fSLois Curfman McInnes 4799b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 480a86d99e1SLois Curfman McInnes 481a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4829b94acceSBarry Smith @*/ 4839b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4849b94acceSBarry Smith { 4853a40ed3dSBarry Smith PetscFunctionBegin; 48677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 48774679c65SBarry Smith PetscValidScalarPointer(gnorm); 48874679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 489d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 49074679c65SBarry Smith } 4919b94acceSBarry Smith *gnorm = snes->norm; 4923a40ed3dSBarry Smith PetscFunctionReturn(0); 4939b94acceSBarry Smith } 49474679c65SBarry Smith 4955615d1e5SSatish Balay #undef __FUNC__ 496d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 4979b94acceSBarry Smith /*@ 4989b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 4999b94acceSBarry Smith attempted by the nonlinear solver. 5009b94acceSBarry Smith 501c7afd0dbSLois Curfman McInnes Not Collective 502c7afd0dbSLois Curfman McInnes 5039b94acceSBarry Smith Input Parameter: 5049b94acceSBarry Smith . snes - SNES context 5059b94acceSBarry Smith 5069b94acceSBarry Smith Output Parameter: 5079b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 5089b94acceSBarry Smith 509c96a6f78SLois Curfman McInnes Notes: 510c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 511c96a6f78SLois Curfman McInnes 51236851e7fSLois Curfman McInnes Level: intermediate 51336851e7fSLois Curfman McInnes 5149b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 5159b94acceSBarry Smith @*/ 5169b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 5179b94acceSBarry Smith { 5183a40ed3dSBarry Smith PetscFunctionBegin; 51977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 52074679c65SBarry Smith PetscValidIntPointer(nfails); 5219b94acceSBarry Smith *nfails = snes->nfailures; 5223a40ed3dSBarry Smith PetscFunctionReturn(0); 5239b94acceSBarry Smith } 524a847f771SSatish Balay 5255615d1e5SSatish Balay #undef __FUNC__ 526d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" 527c96a6f78SLois Curfman McInnes /*@ 528c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 529c96a6f78SLois Curfman McInnes used by the nonlinear solver. 530c96a6f78SLois Curfman McInnes 531c7afd0dbSLois Curfman McInnes Not Collective 532c7afd0dbSLois Curfman McInnes 533c96a6f78SLois Curfman McInnes Input Parameter: 534c96a6f78SLois Curfman McInnes . snes - SNES context 535c96a6f78SLois Curfman McInnes 536c96a6f78SLois Curfman McInnes Output Parameter: 537c96a6f78SLois Curfman McInnes . lits - number of linear iterations 538c96a6f78SLois Curfman McInnes 539c96a6f78SLois Curfman McInnes Notes: 540c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 541c96a6f78SLois Curfman McInnes 54236851e7fSLois Curfman McInnes Level: intermediate 54336851e7fSLois Curfman McInnes 544c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 545c96a6f78SLois Curfman McInnes @*/ 54686bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 547c96a6f78SLois Curfman McInnes { 5483a40ed3dSBarry Smith PetscFunctionBegin; 549c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 550c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 551c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 5523a40ed3dSBarry Smith PetscFunctionReturn(0); 553c96a6f78SLois Curfman McInnes } 554c96a6f78SLois Curfman McInnes 555c96a6f78SLois Curfman McInnes #undef __FUNC__ 556d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES" 5579b94acceSBarry Smith /*@C 5589b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 5599b94acceSBarry Smith 560c7afd0dbSLois Curfman McInnes Not Collective, but if SNES object is parallel, then SLES object is parallel 561c7afd0dbSLois Curfman McInnes 5629b94acceSBarry Smith Input Parameter: 5639b94acceSBarry Smith . snes - the SNES context 5649b94acceSBarry Smith 5659b94acceSBarry Smith Output Parameter: 5669b94acceSBarry Smith . sles - the SLES context 5679b94acceSBarry Smith 5689b94acceSBarry Smith Notes: 5699b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 5709b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5719b94acceSBarry Smith KSP and PC contexts as well. 5729b94acceSBarry Smith 57336851e7fSLois Curfman McInnes Level: beginner 57436851e7fSLois Curfman McInnes 5759b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 5769b94acceSBarry Smith 5779b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 5789b94acceSBarry Smith @*/ 5799b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 5809b94acceSBarry Smith { 5813a40ed3dSBarry Smith PetscFunctionBegin; 58277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5839b94acceSBarry Smith *sles = snes->sles; 5843a40ed3dSBarry Smith PetscFunctionReturn(0); 5859b94acceSBarry Smith } 58682bf6240SBarry Smith 587e24b481bSBarry Smith #undef __FUNC__ 588e24b481bSBarry Smith #define __FUNC__ "SNESPublish_Petsc" 589e24b481bSBarry Smith static int SNESPublish_Petsc(PetscObject object) 590e24b481bSBarry Smith { 591e24b481bSBarry Smith #if defined(HAVE_AMS) 592e24b481bSBarry Smith SNES v = (SNES) object; 593e24b481bSBarry Smith int ierr; 594e24b481bSBarry Smith 595e24b481bSBarry Smith PetscFunctionBegin; 596e24b481bSBarry Smith 597e24b481bSBarry Smith /* if it is already published then return */ 598e24b481bSBarry Smith if (v->amem >=0 ) PetscFunctionReturn(0); 599e24b481bSBarry Smith 6003f1db9ecSBarry Smith ierr = PetscObjectPublishBaseBegin(object);CHKERRQ(ierr); 601e24b481bSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ, 602e24b481bSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 603e24b481bSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ, 604e24b481bSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 605e24b481bSBarry Smith ierr = PetscObjectPublishBaseEnd(object);CHKERRQ(ierr); 606e24b481bSBarry Smith #else 607e24b481bSBarry Smith PetscFunctionBegin; 608e24b481bSBarry Smith #endif 609e24b481bSBarry Smith PetscFunctionReturn(0); 610e24b481bSBarry Smith } 611e24b481bSBarry Smith 6129b94acceSBarry Smith /* -----------------------------------------------------------*/ 6135615d1e5SSatish Balay #undef __FUNC__ 6145615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 6159b94acceSBarry Smith /*@C 6169b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 6179b94acceSBarry Smith 618c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 619c7afd0dbSLois Curfman McInnes 620c7afd0dbSLois Curfman McInnes Input Parameters: 621c7afd0dbSLois Curfman McInnes + comm - MPI communicator 622c7afd0dbSLois Curfman McInnes - type - type of method, either 623c7afd0dbSLois Curfman McInnes SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations) 624c7afd0dbSLois Curfman McInnes or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization) 6259b94acceSBarry Smith 6269b94acceSBarry Smith Output Parameter: 6279b94acceSBarry Smith . outsnes - the new SNES context 6289b94acceSBarry Smith 629c7afd0dbSLois Curfman McInnes Options Database Keys: 630c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 631c7afd0dbSLois Curfman McInnes and no preconditioning matrix 632c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 633c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 634c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 635c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 636c1f60f51SBarry Smith 63736851e7fSLois Curfman McInnes Level: beginner 63836851e7fSLois Curfman McInnes 6399b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 6409b94acceSBarry Smith 64163a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 6429b94acceSBarry Smith @*/ 6434b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 6449b94acceSBarry Smith { 6459b94acceSBarry Smith int ierr; 6469b94acceSBarry Smith SNES snes; 6479b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 64837fcc0dbSBarry Smith 6493a40ed3dSBarry Smith PetscFunctionBegin; 6509b94acceSBarry Smith *outsnes = 0; 651d64ed03dSBarry Smith if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){ 652d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 653d64ed03dSBarry Smith } 6543f1db9ecSBarry Smith PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView); 6559b94acceSBarry Smith PLogObjectCreate(snes); 656e24b481bSBarry Smith snes->bops->publish = SNESPublish_Petsc; 6579b94acceSBarry Smith snes->max_its = 50; 6589750a799SBarry Smith snes->max_funcs = 10000; 6599b94acceSBarry Smith snes->norm = 0.0; 6605a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 661b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 662b18e04deSLois Curfman McInnes snes->ttol = 0.0; 6639b94acceSBarry Smith snes->atol = 1.e-10; 664d64ed03dSBarry Smith } else { 665b4874afaSBarry Smith snes->rtol = 1.e-8; 666b4874afaSBarry Smith snes->ttol = 0.0; 667b4874afaSBarry Smith snes->atol = 1.e-50; 668b4874afaSBarry Smith } 6699b94acceSBarry Smith snes->xtol = 1.e-8; 670b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 6719b94acceSBarry Smith snes->nfuncs = 0; 6729b94acceSBarry Smith snes->nfailures = 0; 6737a00f4a9SLois Curfman McInnes snes->linear_its = 0; 674639f9d9dSBarry Smith snes->numbermonitors = 0; 6759b94acceSBarry Smith snes->data = 0; 6769b94acceSBarry Smith snes->view = 0; 6779b94acceSBarry Smith snes->computeumfunction = 0; 6789b94acceSBarry Smith snes->umfunP = 0; 6799b94acceSBarry Smith snes->fc = 0; 6809b94acceSBarry Smith snes->deltatol = 1.e-12; 6819b94acceSBarry Smith snes->fmin = -1.e30; 6829b94acceSBarry Smith snes->method_class = type; 6839b94acceSBarry Smith snes->set_method_called = 0; 68482bf6240SBarry Smith snes->setupcalled = 0; 6859b94acceSBarry Smith snes->ksp_ewconv = 0; 6866f24a144SLois Curfman McInnes snes->xmonitor = 0; 6876f24a144SLois Curfman McInnes snes->vwork = 0; 6886f24a144SLois Curfman McInnes snes->nwork = 0; 689758f92a0SBarry Smith snes->conv_hist_len = 0; 690758f92a0SBarry Smith snes->conv_hist_max = 0; 691758f92a0SBarry Smith snes->conv_hist = PETSC_NULL; 692758f92a0SBarry Smith snes->conv_hist_its = PETSC_NULL; 693758f92a0SBarry Smith snes->conv_hist_reset = PETSC_TRUE; 6949b94acceSBarry Smith 6959b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 6960452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx); CHKPTRQ(kctx); 697eed86810SBarry Smith PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 6989b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 6999b94acceSBarry Smith kctx->version = 2; 7009b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 7019b94acceSBarry Smith this was too large for some test cases */ 7029b94acceSBarry Smith kctx->rtol_last = 0; 7039b94acceSBarry Smith kctx->rtol_max = .9; 7049b94acceSBarry Smith kctx->gamma = 1.0; 7059b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 7069b94acceSBarry Smith kctx->alpha = kctx->alpha2; 7079b94acceSBarry Smith kctx->threshold = .1; 7089b94acceSBarry Smith kctx->lresid_last = 0; 7099b94acceSBarry Smith kctx->norm_last = 0; 7109b94acceSBarry Smith 7119b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles); CHKERRQ(ierr); 7129b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 7135334005bSBarry Smith 7149b94acceSBarry Smith *outsnes = snes; 715e24b481bSBarry Smith PetscPublishAll(snes); 7163a40ed3dSBarry Smith PetscFunctionReturn(0); 7179b94acceSBarry Smith } 7189b94acceSBarry Smith 7199b94acceSBarry Smith /* --------------------------------------------------------------- */ 7205615d1e5SSatish Balay #undef __FUNC__ 721d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction" 7229b94acceSBarry Smith /*@C 7239b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 7249b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 7259b94acceSBarry Smith equations. 7269b94acceSBarry Smith 727fee21e36SBarry Smith Collective on SNES 728fee21e36SBarry Smith 729c7afd0dbSLois Curfman McInnes Input Parameters: 730c7afd0dbSLois Curfman McInnes + snes - the SNES context 731c7afd0dbSLois Curfman McInnes . func - function evaluation routine 732c7afd0dbSLois Curfman McInnes . r - vector to store function value 733c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 734c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 7359b94acceSBarry Smith 736c7afd0dbSLois Curfman McInnes Calling sequence of func: 7378d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 738c7afd0dbSLois Curfman McInnes 739313e4042SLois Curfman McInnes . f - function vector 740c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 7419b94acceSBarry Smith 7429b94acceSBarry Smith Notes: 7439b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 7449b94acceSBarry Smith $ f'(x) x = -f(x), 745c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 7469b94acceSBarry Smith 7479b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 7489b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 7499b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 7509b94acceSBarry Smith 75136851e7fSLois Curfman McInnes Level: beginner 75236851e7fSLois Curfman McInnes 7539b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 7549b94acceSBarry Smith 755a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 7569b94acceSBarry Smith @*/ 7575334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 7589b94acceSBarry Smith { 7593a40ed3dSBarry Smith PetscFunctionBegin; 76077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 761a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 762a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 763a8c6a408SBarry Smith } 7649b94acceSBarry Smith snes->computefunction = func; 7659b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7669b94acceSBarry Smith snes->funP = ctx; 7673a40ed3dSBarry Smith PetscFunctionReturn(0); 7689b94acceSBarry Smith } 7699b94acceSBarry Smith 7705615d1e5SSatish Balay #undef __FUNC__ 7715615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 7729b94acceSBarry Smith /*@ 77336851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7749b94acceSBarry Smith SNESSetFunction(). 7759b94acceSBarry Smith 776c7afd0dbSLois Curfman McInnes Collective on SNES 777c7afd0dbSLois Curfman McInnes 7789b94acceSBarry Smith Input Parameters: 779c7afd0dbSLois Curfman McInnes + snes - the SNES context 780c7afd0dbSLois Curfman McInnes - x - input vector 7819b94acceSBarry Smith 7829b94acceSBarry Smith Output Parameter: 7833638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7849b94acceSBarry Smith 7851bffabb2SLois Curfman McInnes Notes: 7869b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 7879b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 7889b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 7899b94acceSBarry Smith 79036851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 79136851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 79236851e7fSLois Curfman McInnes themselves. 79336851e7fSLois Curfman McInnes 79436851e7fSLois Curfman McInnes Level: developer 79536851e7fSLois Curfman McInnes 7969b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 7979b94acceSBarry Smith 798a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 7999b94acceSBarry Smith @*/ 8009b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 8019b94acceSBarry Smith { 8029b94acceSBarry Smith int ierr; 8039b94acceSBarry Smith 8043a40ed3dSBarry Smith PetscFunctionBegin; 805d4bb536fSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 806a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 807d4bb536fSBarry Smith } 8089b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 809d64ed03dSBarry Smith PetscStackPush("SNES user function"); 8109b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 811d64ed03dSBarry Smith PetscStackPop; 812ae3c334cSLois Curfman McInnes snes->nfuncs++; 8139b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 8143a40ed3dSBarry Smith PetscFunctionReturn(0); 8159b94acceSBarry Smith } 8169b94acceSBarry Smith 8175615d1e5SSatish Balay #undef __FUNC__ 818d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" 8199b94acceSBarry Smith /*@C 8209b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 8219b94acceSBarry Smith unconstrained minimization. 8229b94acceSBarry Smith 823fee21e36SBarry Smith Collective on SNES 824fee21e36SBarry Smith 825c7afd0dbSLois Curfman McInnes Input Parameters: 826c7afd0dbSLois Curfman McInnes + snes - the SNES context 827c7afd0dbSLois Curfman McInnes . func - function evaluation routine 828c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 829c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 8309b94acceSBarry Smith 831c7afd0dbSLois Curfman McInnes Calling sequence of func: 8328d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,double *f,void *ctx); 833c7afd0dbSLois Curfman McInnes 834c7afd0dbSLois Curfman McInnes + x - input vector 8359b94acceSBarry Smith . f - function 836c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined function context 8379b94acceSBarry Smith 83836851e7fSLois Curfman McInnes Level: beginner 83936851e7fSLois Curfman McInnes 8409b94acceSBarry Smith Notes: 8419b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 8429b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 8439b94acceSBarry Smith SNESSetFunction(). 8449b94acceSBarry Smith 8459b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 8469b94acceSBarry Smith 847a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 848a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 8499b94acceSBarry Smith @*/ 8509b94acceSBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*), 8519b94acceSBarry Smith void *ctx) 8529b94acceSBarry Smith { 8533a40ed3dSBarry Smith PetscFunctionBegin; 85477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 855a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 856a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 857a8c6a408SBarry Smith } 8589b94acceSBarry Smith snes->computeumfunction = func; 8599b94acceSBarry Smith snes->umfunP = ctx; 8603a40ed3dSBarry Smith PetscFunctionReturn(0); 8619b94acceSBarry Smith } 8629b94acceSBarry Smith 8635615d1e5SSatish Balay #undef __FUNC__ 8645615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 8659b94acceSBarry Smith /*@ 8669b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 8679b94acceSBarry Smith set with SNESSetMinimizationFunction(). 8689b94acceSBarry Smith 869c7afd0dbSLois Curfman McInnes Collective on SNES 870c7afd0dbSLois Curfman McInnes 8719b94acceSBarry Smith Input Parameters: 872c7afd0dbSLois Curfman McInnes + snes - the SNES context 873c7afd0dbSLois Curfman McInnes - x - input vector 8749b94acceSBarry Smith 8759b94acceSBarry Smith Output Parameter: 8769b94acceSBarry Smith . y - function value 8779b94acceSBarry Smith 8789b94acceSBarry Smith Notes: 8799b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 8809b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 8819b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 882a86d99e1SLois Curfman McInnes 88336851e7fSLois Curfman McInnes SNESComputeMinimizationFunction() is typically used within minimization 88436851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 88536851e7fSLois Curfman McInnes themselves. 88636851e7fSLois Curfman McInnes 88736851e7fSLois Curfman McInnes Level: developer 88836851e7fSLois Curfman McInnes 889a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 890a86d99e1SLois Curfman McInnes 891a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 892a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 8939b94acceSBarry Smith @*/ 8949b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 8959b94acceSBarry Smith { 8969b94acceSBarry Smith int ierr; 8973a40ed3dSBarry Smith 8983a40ed3dSBarry Smith PetscFunctionBegin; 899a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 900a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 901a8c6a408SBarry Smith } 9029b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 903d64ed03dSBarry Smith PetscStackPush("SNES user minimzation function"); 9049b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP); CHKERRQ(ierr); 905d64ed03dSBarry Smith PetscStackPop; 906ae3c334cSLois Curfman McInnes snes->nfuncs++; 9079b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 9083a40ed3dSBarry Smith PetscFunctionReturn(0); 9099b94acceSBarry Smith } 9109b94acceSBarry Smith 9115615d1e5SSatish Balay #undef __FUNC__ 912d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient" 9139b94acceSBarry Smith /*@C 9149b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 9159b94acceSBarry Smith vector for use by the SNES routines. 9169b94acceSBarry Smith 917c7afd0dbSLois Curfman McInnes Collective on SNES 918c7afd0dbSLois Curfman McInnes 9199b94acceSBarry Smith Input Parameters: 920c7afd0dbSLois Curfman McInnes + snes - the SNES context 9219b94acceSBarry Smith . func - function evaluation routine 922044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 923044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 924c7afd0dbSLois Curfman McInnes - r - vector to store gradient value 925fee21e36SBarry Smith 9269b94acceSBarry Smith Calling sequence of func: 9278d76a1e5SLois Curfman McInnes $ func (SNES, Vec x, Vec g, void *ctx); 9289b94acceSBarry Smith 929c7afd0dbSLois Curfman McInnes + x - input vector 9309b94acceSBarry Smith . g - gradient vector 931c7afd0dbSLois Curfman McInnes - ctx - optional user-defined gradient context 9329b94acceSBarry Smith 9339b94acceSBarry Smith Notes: 9349b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 9359b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 9369b94acceSBarry Smith SNESSetFunction(). 9379b94acceSBarry Smith 93836851e7fSLois Curfman McInnes Level: beginner 93936851e7fSLois Curfman McInnes 9409b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 9419b94acceSBarry Smith 942a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 943a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 9449b94acceSBarry Smith @*/ 94574679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 9469b94acceSBarry Smith { 9473a40ed3dSBarry Smith PetscFunctionBegin; 94877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 949a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 950a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 951a8c6a408SBarry Smith } 9529b94acceSBarry Smith snes->computefunction = func; 9539b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 9549b94acceSBarry Smith snes->funP = ctx; 9553a40ed3dSBarry Smith PetscFunctionReturn(0); 9569b94acceSBarry Smith } 9579b94acceSBarry Smith 9585615d1e5SSatish Balay #undef __FUNC__ 9595615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 9609b94acceSBarry Smith /*@ 961a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 962a86d99e1SLois Curfman McInnes SNESSetGradient(). 9639b94acceSBarry Smith 964c7afd0dbSLois Curfman McInnes Collective on SNES 965c7afd0dbSLois Curfman McInnes 9669b94acceSBarry Smith Input Parameters: 967c7afd0dbSLois Curfman McInnes + snes - the SNES context 968c7afd0dbSLois Curfman McInnes - x - input vector 9699b94acceSBarry Smith 9709b94acceSBarry Smith Output Parameter: 9719b94acceSBarry Smith . y - gradient vector 9729b94acceSBarry Smith 9739b94acceSBarry Smith Notes: 9749b94acceSBarry Smith SNESComputeGradient() is valid only for 9759b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 9769b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 977a86d99e1SLois Curfman McInnes 97836851e7fSLois Curfman McInnes SNESComputeGradient() is typically used within minimization 97936851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 98036851e7fSLois Curfman McInnes themselves. 98136851e7fSLois Curfman McInnes 98236851e7fSLois Curfman McInnes Level: developer 98336851e7fSLois Curfman McInnes 984a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 985a86d99e1SLois Curfman McInnes 986a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 987a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 9889b94acceSBarry Smith @*/ 9899b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 9909b94acceSBarry Smith { 9919b94acceSBarry Smith int ierr; 9923a40ed3dSBarry Smith 9933a40ed3dSBarry Smith PetscFunctionBegin; 9943a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 995a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 9963a40ed3dSBarry Smith } 9979b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 998d64ed03dSBarry Smith PetscStackPush("SNES user gradient function"); 9999b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP); CHKERRQ(ierr); 1000d64ed03dSBarry Smith PetscStackPop; 10019b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 10023a40ed3dSBarry Smith PetscFunctionReturn(0); 10039b94acceSBarry Smith } 10049b94acceSBarry Smith 10055615d1e5SSatish Balay #undef __FUNC__ 10065615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 100762fef451SLois Curfman McInnes /*@ 100862fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 100962fef451SLois Curfman McInnes set with SNESSetJacobian(). 101062fef451SLois Curfman McInnes 1011c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1012c7afd0dbSLois Curfman McInnes 101362fef451SLois Curfman McInnes Input Parameters: 1014c7afd0dbSLois Curfman McInnes + snes - the SNES context 1015c7afd0dbSLois Curfman McInnes - x - input vector 101662fef451SLois Curfman McInnes 101762fef451SLois Curfman McInnes Output Parameters: 1018c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 101962fef451SLois Curfman McInnes . B - optional preconditioning matrix 1020c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 1021fee21e36SBarry Smith 102262fef451SLois Curfman McInnes Notes: 102362fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 102462fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 102562fef451SLois Curfman McInnes 1026dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 1027dc5a77f8SLois Curfman McInnes flag parameter. 102862fef451SLois Curfman McInnes 102962fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 103062fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 1031005c665bSBarry Smith methods is SNESComputeHessian(). 103262fef451SLois Curfman McInnes 103336851e7fSLois Curfman McInnes SNESComputeJacobian() is typically used within nonlinear solver 103436851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 103536851e7fSLois Curfman McInnes themselves. 103636851e7fSLois Curfman McInnes 103736851e7fSLois Curfman McInnes Level: developer 103836851e7fSLois Curfman McInnes 103962fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 104062fef451SLois Curfman McInnes 104162fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 104262fef451SLois Curfman McInnes @*/ 10439b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 10449b94acceSBarry Smith { 10459b94acceSBarry Smith int ierr; 10463a40ed3dSBarry Smith 10473a40ed3dSBarry Smith PetscFunctionBegin; 10483a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1049a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 10503a40ed3dSBarry Smith } 10513a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 10529b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 1053c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 1054d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 10559b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP); CHKERRQ(ierr); 1056d64ed03dSBarry Smith PetscStackPop; 10579b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 10586d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 105977c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 106077c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 10613a40ed3dSBarry Smith PetscFunctionReturn(0); 10629b94acceSBarry Smith } 10639b94acceSBarry Smith 10645615d1e5SSatish Balay #undef __FUNC__ 10655615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 106662fef451SLois Curfman McInnes /*@ 106762fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 106862fef451SLois Curfman McInnes set with SNESSetHessian(). 106962fef451SLois Curfman McInnes 1070c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1071c7afd0dbSLois Curfman McInnes 107262fef451SLois Curfman McInnes Input Parameters: 1073c7afd0dbSLois Curfman McInnes + snes - the SNES context 1074c7afd0dbSLois Curfman McInnes - x - input vector 107562fef451SLois Curfman McInnes 107662fef451SLois Curfman McInnes Output Parameters: 1077c7afd0dbSLois Curfman McInnes + A - Hessian matrix 107862fef451SLois Curfman McInnes . B - optional preconditioning matrix 1079c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 1080fee21e36SBarry Smith 108162fef451SLois Curfman McInnes Notes: 108262fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 108362fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 108462fef451SLois Curfman McInnes 1085dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 1086dc5a77f8SLois Curfman McInnes flag parameter. 108762fef451SLois Curfman McInnes 108862fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 108962fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 109062fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 109162fef451SLois Curfman McInnes 109236851e7fSLois Curfman McInnes SNESComputeHessian() is typically used within minimization 109336851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 109436851e7fSLois Curfman McInnes themselves. 109536851e7fSLois Curfman McInnes 109636851e7fSLois Curfman McInnes Level: developer 109736851e7fSLois Curfman McInnes 109862fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 109962fef451SLois Curfman McInnes 1100a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 1101a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 110262fef451SLois Curfman McInnes @*/ 110362fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 11049b94acceSBarry Smith { 11059b94acceSBarry Smith int ierr; 11063a40ed3dSBarry Smith 11073a40ed3dSBarry Smith PetscFunctionBegin; 11083a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1109a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 11103a40ed3dSBarry Smith } 11113a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 111262fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 11130f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 1114d64ed03dSBarry Smith PetscStackPush("SNES user Hessian function"); 111562fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP); CHKERRQ(ierr); 1116d64ed03dSBarry Smith PetscStackPop; 111762fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 11180f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 111977c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 112077c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 11213a40ed3dSBarry Smith PetscFunctionReturn(0); 11229b94acceSBarry Smith } 11239b94acceSBarry Smith 11245615d1e5SSatish Balay #undef __FUNC__ 1125d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian" 11269b94acceSBarry Smith /*@C 11279b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 1128044dda88SLois Curfman McInnes location to store the matrix. 11299b94acceSBarry Smith 1130c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1131c7afd0dbSLois Curfman McInnes 11329b94acceSBarry Smith Input Parameters: 1133c7afd0dbSLois Curfman McInnes + snes - the SNES context 11349b94acceSBarry Smith . A - Jacobian matrix 11359b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 11369b94acceSBarry Smith . func - Jacobian evaluation routine 1137c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 11382cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 11399b94acceSBarry Smith 11409b94acceSBarry Smith Calling sequence of func: 11418d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 11429b94acceSBarry Smith 1143c7afd0dbSLois Curfman McInnes + x - input vector 11449b94acceSBarry Smith . A - Jacobian matrix 11459b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1146ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1147ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 1148c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 11499b94acceSBarry Smith 11509b94acceSBarry Smith Notes: 1151dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 11522cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1153ac21db08SLois Curfman McInnes 1154ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 11559b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 11569b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 11579b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 11589b94acceSBarry Smith throughout the global iterations. 11599b94acceSBarry Smith 116036851e7fSLois Curfman McInnes Level: beginner 116136851e7fSLois Curfman McInnes 11629b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 11639b94acceSBarry Smith 1164ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 11659b94acceSBarry Smith @*/ 11669b94acceSBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 11679b94acceSBarry Smith MatStructure*,void*),void *ctx) 11689b94acceSBarry Smith { 11693a40ed3dSBarry Smith PetscFunctionBegin; 117077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1171a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1172a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1173a8c6a408SBarry Smith } 11749b94acceSBarry Smith snes->computejacobian = func; 11759b94acceSBarry Smith snes->jacP = ctx; 11769b94acceSBarry Smith snes->jacobian = A; 11779b94acceSBarry Smith snes->jacobian_pre = B; 11783a40ed3dSBarry Smith PetscFunctionReturn(0); 11799b94acceSBarry Smith } 118062fef451SLois Curfman McInnes 11815615d1e5SSatish Balay #undef __FUNC__ 1182d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian" 1183b4fd4287SBarry Smith /*@ 1184b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1185b4fd4287SBarry Smith provided context for evaluating the Jacobian. 1186b4fd4287SBarry Smith 1187c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 1188c7afd0dbSLois Curfman McInnes 1189b4fd4287SBarry Smith Input Parameter: 1190b4fd4287SBarry Smith . snes - the nonlinear solver context 1191b4fd4287SBarry Smith 1192b4fd4287SBarry Smith Output Parameters: 1193c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 1194b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 1195c7afd0dbSLois Curfman McInnes - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1196fee21e36SBarry Smith 119736851e7fSLois Curfman McInnes Level: advanced 119836851e7fSLois Curfman McInnes 1199b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 1200b4fd4287SBarry Smith @*/ 1201b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1202b4fd4287SBarry Smith { 12033a40ed3dSBarry Smith PetscFunctionBegin; 12043a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1205a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 12063a40ed3dSBarry Smith } 1207b4fd4287SBarry Smith if (A) *A = snes->jacobian; 1208b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 1209b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 12103a40ed3dSBarry Smith PetscFunctionReturn(0); 1211b4fd4287SBarry Smith } 1212b4fd4287SBarry Smith 12135615d1e5SSatish Balay #undef __FUNC__ 1214d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian" 12159b94acceSBarry Smith /*@C 12169b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 1217044dda88SLois Curfman McInnes location to store the matrix. 12189b94acceSBarry Smith 1219c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1220c7afd0dbSLois Curfman McInnes 12219b94acceSBarry Smith Input Parameters: 1222c7afd0dbSLois Curfman McInnes + snes - the SNES context 12239b94acceSBarry Smith . A - Hessian matrix 12249b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 12259b94acceSBarry Smith . func - Jacobian evaluation routine 1226c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 1227313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 12289b94acceSBarry Smith 12299b94acceSBarry Smith Calling sequence of func: 12308d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 12319b94acceSBarry Smith 1232c7afd0dbSLois Curfman McInnes + x - input vector 12339b94acceSBarry Smith . A - Hessian matrix 12349b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1235ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1236ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 1237c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Hessian context 12389b94acceSBarry Smith 12399b94acceSBarry Smith Notes: 1240dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 12412cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1242ac21db08SLois Curfman McInnes 12439b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 12449b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 12459b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 12469b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 12479b94acceSBarry Smith throughout the global iterations. 12489b94acceSBarry Smith 124936851e7fSLois Curfman McInnes Level: beginner 125036851e7fSLois Curfman McInnes 12519b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 12529b94acceSBarry Smith 1253ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 12549b94acceSBarry Smith @*/ 12559b94acceSBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*, 12569b94acceSBarry Smith MatStructure*,void*),void *ctx) 12579b94acceSBarry Smith { 12583a40ed3dSBarry Smith PetscFunctionBegin; 125977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1260d4bb536fSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1261a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1262d4bb536fSBarry Smith } 12639b94acceSBarry Smith snes->computejacobian = func; 12649b94acceSBarry Smith snes->jacP = ctx; 12659b94acceSBarry Smith snes->jacobian = A; 12669b94acceSBarry Smith snes->jacobian_pre = B; 12673a40ed3dSBarry Smith PetscFunctionReturn(0); 12689b94acceSBarry Smith } 12699b94acceSBarry Smith 12705615d1e5SSatish Balay #undef __FUNC__ 1271d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian" 127262fef451SLois Curfman McInnes /*@ 127362fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 127462fef451SLois Curfman McInnes provided context for evaluating the Hessian. 127562fef451SLois Curfman McInnes 1276c7afd0dbSLois Curfman McInnes Not Collective, but Mat object is parallel if SNES object is parallel 1277c7afd0dbSLois Curfman McInnes 127862fef451SLois Curfman McInnes Input Parameter: 127962fef451SLois Curfman McInnes . snes - the nonlinear solver context 128062fef451SLois Curfman McInnes 128162fef451SLois Curfman McInnes Output Parameters: 1282c7afd0dbSLois Curfman McInnes + A - location to stash Hessian matrix (or PETSC_NULL) 128362fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 1284c7afd0dbSLois Curfman McInnes - ctx - location to stash Hessian ctx (or PETSC_NULL) 1285fee21e36SBarry Smith 128636851e7fSLois Curfman McInnes Level: advanced 128736851e7fSLois Curfman McInnes 128862fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 1289c7afd0dbSLois Curfman McInnes 1290c7afd0dbSLois Curfman McInnes .keywords: SNES, get, Hessian 129162fef451SLois Curfman McInnes @*/ 129262fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 129362fef451SLois Curfman McInnes { 12943a40ed3dSBarry Smith PetscFunctionBegin; 12953a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){ 1296a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 12973a40ed3dSBarry Smith } 129862fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 129962fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 130062fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 13013a40ed3dSBarry Smith PetscFunctionReturn(0); 130262fef451SLois Curfman McInnes } 130362fef451SLois Curfman McInnes 13049b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 13059b94acceSBarry Smith 13065615d1e5SSatish Balay #undef __FUNC__ 13075615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 13089b94acceSBarry Smith /*@ 13099b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1310272ac6f2SLois Curfman McInnes of a nonlinear solver. 13119b94acceSBarry Smith 1312fee21e36SBarry Smith Collective on SNES 1313fee21e36SBarry Smith 1314c7afd0dbSLois Curfman McInnes Input Parameters: 1315c7afd0dbSLois Curfman McInnes + snes - the SNES context 1316c7afd0dbSLois Curfman McInnes - x - the solution vector 1317c7afd0dbSLois Curfman McInnes 1318272ac6f2SLois Curfman McInnes Notes: 1319272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1320272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1321272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1322272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1323272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1324272ac6f2SLois Curfman McInnes 132536851e7fSLois Curfman McInnes Level: advanced 132636851e7fSLois Curfman McInnes 13279b94acceSBarry Smith .keywords: SNES, nonlinear, setup 13289b94acceSBarry Smith 13299b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 13309b94acceSBarry Smith @*/ 13318ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 13329b94acceSBarry Smith { 1333272ac6f2SLois Curfman McInnes int ierr, flg; 13343a40ed3dSBarry Smith 13353a40ed3dSBarry Smith PetscFunctionBegin; 133677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 133777c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 13388ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 13399b94acceSBarry Smith 1340c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg); CHKERRQ(ierr); 1341c1f60f51SBarry Smith /* 1342c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1343dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1344c1f60f51SBarry Smith */ 1345c1f60f51SBarry Smith if (flg) { 1346c1f60f51SBarry Smith Mat J; 13475a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1348c1f60f51SBarry Smith PLogObjectParent(snes,J); 1349c1f60f51SBarry Smith snes->mfshell = J; 1350c1f60f51SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1351c1f60f51SBarry Smith snes->jacobian = J; 1352c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1353d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1354c1f60f51SBarry Smith snes->jacobian = J; 1355c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1356d4bb536fSBarry Smith } else { 1357a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option"); 1358d4bb536fSBarry Smith } 13595a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1360c1f60f51SBarry Smith } 1361272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg); CHKERRQ(ierr); 1362c1f60f51SBarry Smith /* 1363dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1364c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1365c1f60f51SBarry Smith */ 136631615d04SLois Curfman McInnes if (flg) { 1367272ac6f2SLois Curfman McInnes Mat J; 13685a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1369272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1370272ac6f2SLois Curfman McInnes snes->mfshell = J; 137183e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 137283e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 137394a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1374d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 137583e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP); CHKERRQ(ierr); 137694a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1377d4bb536fSBarry Smith } else { 1378a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option"); 1379d4bb536fSBarry Smith } 13805a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1381272ac6f2SLois Curfman McInnes } 13829b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1383a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1384a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 13855a655dc6SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1386a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 1387a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector"); 1388a8c6a408SBarry Smith } 1389a703fe33SLois Curfman McInnes 1390a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 139182bf6240SBarry Smith if (snes->ksp_ewconv && PetscStrcmp(snes->type_name,SNES_EQ_TR)) { 1392a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1393a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles); CHKERRQ(ierr); 1394a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp); CHKERRQ(ierr); 13955a655dc6SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,(void *)snes); CHKERRQ(ierr); 1396a703fe33SLois Curfman McInnes } 13979b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1398a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1399a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1400a8c6a408SBarry Smith if (!snes->computeumfunction) { 1401a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first"); 1402a8c6a408SBarry Smith } 14035a655dc6SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian()"); 1404d4bb536fSBarry Smith } else { 1405a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class"); 1406d4bb536fSBarry Smith } 1407a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes); CHKERRQ(ierr);} 140882bf6240SBarry Smith snes->setupcalled = 1; 14093a40ed3dSBarry Smith PetscFunctionReturn(0); 14109b94acceSBarry Smith } 14119b94acceSBarry Smith 14125615d1e5SSatish Balay #undef __FUNC__ 1413d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy" 14149b94acceSBarry Smith /*@C 14159b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 14169b94acceSBarry Smith with SNESCreate(). 14179b94acceSBarry Smith 1418c7afd0dbSLois Curfman McInnes Collective on SNES 1419c7afd0dbSLois Curfman McInnes 14209b94acceSBarry Smith Input Parameter: 14219b94acceSBarry Smith . snes - the SNES context 14229b94acceSBarry Smith 142336851e7fSLois Curfman McInnes Level: beginner 142436851e7fSLois Curfman McInnes 14259b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 14269b94acceSBarry Smith 142763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 14289b94acceSBarry Smith @*/ 14299b94acceSBarry Smith int SNESDestroy(SNES snes) 14309b94acceSBarry Smith { 14319b94acceSBarry Smith int ierr; 14323a40ed3dSBarry Smith 14333a40ed3dSBarry Smith PetscFunctionBegin; 143477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 14353a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1436d4bb536fSBarry Smith 1437e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes); CHKERRQ(ierr);} 14380452661fSBarry Smith if (snes->kspconvctx) PetscFree(snes->kspconvctx); 1439522c5e43SBarry Smith if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);} 14409b94acceSBarry Smith ierr = SLESDestroy(snes->sles); CHKERRQ(ierr); 1441522c5e43SBarry Smith if (snes->xmonitor) {ierr = SNESLGMonitorDestroy(snes->xmonitor);CHKERRQ(ierr);} 1442522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 14439b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 14440452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 14453a40ed3dSBarry Smith PetscFunctionReturn(0); 14469b94acceSBarry Smith } 14479b94acceSBarry Smith 14489b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 14499b94acceSBarry Smith 14505615d1e5SSatish Balay #undef __FUNC__ 14515615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 14529b94acceSBarry Smith /*@ 1453d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 14549b94acceSBarry Smith 1455c7afd0dbSLois Curfman McInnes Collective on SNES 1456c7afd0dbSLois Curfman McInnes 14579b94acceSBarry Smith Input Parameters: 1458c7afd0dbSLois Curfman McInnes + snes - the SNES context 145933174efeSLois Curfman McInnes . atol - absolute convergence tolerance 146033174efeSLois Curfman McInnes . rtol - relative convergence tolerance 146133174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 146233174efeSLois Curfman McInnes of the change in the solution between steps 146333174efeSLois Curfman McInnes . maxit - maximum number of iterations 1464c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1465fee21e36SBarry Smith 146633174efeSLois Curfman McInnes Options Database Keys: 1467c7afd0dbSLois Curfman McInnes + -snes_atol <atol> - Sets atol 1468c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1469c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1470c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1471c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 14729b94acceSBarry Smith 1473d7a720efSLois Curfman McInnes Notes: 14749b94acceSBarry Smith The default maximum number of iterations is 50. 14759b94acceSBarry Smith The default maximum number of function evaluations is 1000. 14769b94acceSBarry Smith 147736851e7fSLois Curfman McInnes Level: intermediate 147836851e7fSLois Curfman McInnes 147933174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 14809b94acceSBarry Smith 1481d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 14829b94acceSBarry Smith @*/ 1483d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 14849b94acceSBarry Smith { 14853a40ed3dSBarry Smith PetscFunctionBegin; 148677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1487d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1488d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1489d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1490d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1491d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 14923a40ed3dSBarry Smith PetscFunctionReturn(0); 14939b94acceSBarry Smith } 14949b94acceSBarry Smith 14955615d1e5SSatish Balay #undef __FUNC__ 14965615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 14979b94acceSBarry Smith /*@ 149833174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 149933174efeSLois Curfman McInnes 1500c7afd0dbSLois Curfman McInnes Not Collective 1501c7afd0dbSLois Curfman McInnes 150233174efeSLois Curfman McInnes Input Parameters: 1503c7afd0dbSLois Curfman McInnes + snes - the SNES context 150433174efeSLois Curfman McInnes . atol - absolute convergence tolerance 150533174efeSLois Curfman McInnes . rtol - relative convergence tolerance 150633174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 150733174efeSLois Curfman McInnes of the change in the solution between steps 150833174efeSLois Curfman McInnes . maxit - maximum number of iterations 1509c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1510fee21e36SBarry Smith 151133174efeSLois Curfman McInnes Notes: 151233174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 151333174efeSLois Curfman McInnes 151436851e7fSLois Curfman McInnes Level: intermediate 151536851e7fSLois Curfman McInnes 151633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 151733174efeSLois Curfman McInnes 151833174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 151933174efeSLois Curfman McInnes @*/ 152033174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 152133174efeSLois Curfman McInnes { 15223a40ed3dSBarry Smith PetscFunctionBegin; 152333174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 152433174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 152533174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 152633174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 152733174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 152833174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 15293a40ed3dSBarry Smith PetscFunctionReturn(0); 153033174efeSLois Curfman McInnes } 153133174efeSLois Curfman McInnes 15325615d1e5SSatish Balay #undef __FUNC__ 15335615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 153433174efeSLois Curfman McInnes /*@ 15359b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 15369b94acceSBarry Smith 1537fee21e36SBarry Smith Collective on SNES 1538fee21e36SBarry Smith 1539c7afd0dbSLois Curfman McInnes Input Parameters: 1540c7afd0dbSLois Curfman McInnes + snes - the SNES context 1541c7afd0dbSLois Curfman McInnes - tol - tolerance 1542c7afd0dbSLois Curfman McInnes 15439b94acceSBarry Smith Options Database Key: 1544c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 15459b94acceSBarry Smith 154636851e7fSLois Curfman McInnes Level: intermediate 154736851e7fSLois Curfman McInnes 15489b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 15499b94acceSBarry Smith 1550d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 15519b94acceSBarry Smith @*/ 15529b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 15539b94acceSBarry Smith { 15543a40ed3dSBarry Smith PetscFunctionBegin; 155577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 15569b94acceSBarry Smith snes->deltatol = tol; 15573a40ed3dSBarry Smith PetscFunctionReturn(0); 15589b94acceSBarry Smith } 15599b94acceSBarry Smith 15605615d1e5SSatish Balay #undef __FUNC__ 15615615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 15629b94acceSBarry Smith /*@ 156377c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 15649b94acceSBarry Smith for unconstrained minimization solvers. 15659b94acceSBarry Smith 1566fee21e36SBarry Smith Collective on SNES 1567fee21e36SBarry Smith 1568c7afd0dbSLois Curfman McInnes Input Parameters: 1569c7afd0dbSLois Curfman McInnes + snes - the SNES context 1570c7afd0dbSLois Curfman McInnes - ftol - minimum function tolerance 1571c7afd0dbSLois Curfman McInnes 15729b94acceSBarry Smith Options Database Key: 1573c7afd0dbSLois Curfman McInnes . -snes_fmin <ftol> - Sets ftol 15749b94acceSBarry Smith 15759b94acceSBarry Smith Note: 157677c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 15779b94acceSBarry Smith methods only. 15789b94acceSBarry Smith 157936851e7fSLois Curfman McInnes Level: intermediate 158036851e7fSLois Curfman McInnes 15819b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 15829b94acceSBarry Smith 1583d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 15849b94acceSBarry Smith @*/ 158577c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 15869b94acceSBarry Smith { 15873a40ed3dSBarry Smith PetscFunctionBegin; 158877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 15899b94acceSBarry Smith snes->fmin = ftol; 15903a40ed3dSBarry Smith PetscFunctionReturn(0); 15919b94acceSBarry Smith } 15929b94acceSBarry Smith 1593ce1608b8SBarry Smith #undef __FUNC__ 1594ce1608b8SBarry Smith #define __FUNC__ "SNESLGMonitor" 1595ce1608b8SBarry Smith int SNESLGMonitor(SNES snes,int it,double norm,void *ctx) 1596ce1608b8SBarry Smith { 1597ce1608b8SBarry Smith int ierr; 1598ce1608b8SBarry Smith 1599ce1608b8SBarry Smith PetscFunctionBegin; 1600ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1601ce1608b8SBarry Smith PetscFunctionReturn(0); 1602ce1608b8SBarry Smith } 1603ce1608b8SBarry Smith 16049b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 16059b94acceSBarry Smith 16065615d1e5SSatish Balay #undef __FUNC__ 1607d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor" 16089b94acceSBarry Smith /*@C 16095cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 16109b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 16119b94acceSBarry Smith progress. 16129b94acceSBarry Smith 1613fee21e36SBarry Smith Collective on SNES 1614fee21e36SBarry Smith 1615c7afd0dbSLois Curfman McInnes Input Parameters: 1616c7afd0dbSLois Curfman McInnes + snes - the SNES context 1617c7afd0dbSLois Curfman McInnes . func - monitoring routine 1618c7afd0dbSLois Curfman McInnes - mctx - [optional] user-defined context for private data for the 1619c7afd0dbSLois Curfman McInnes monitor routine (may be PETSC_NULL) 16209b94acceSBarry Smith 1621c7afd0dbSLois Curfman McInnes Calling sequence of func: 162240a0c1c6SLois Curfman McInnes $ int func(SNES snes,int its, double norm,void *mctx) 1623c7afd0dbSLois Curfman McInnes 1624c7afd0dbSLois Curfman McInnes + snes - the SNES context 1625c7afd0dbSLois Curfman McInnes . its - iteration number 1626c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 1627c7afd0dbSLois Curfman McInnes for SNES_NONLINEAR_EQUATIONS methods 162840a0c1c6SLois Curfman McInnes . norm - 2-norm gradient value (may be estimated) 1629c7afd0dbSLois Curfman McInnes for SNES_UNCONSTRAINED_MINIMIZATION methods 163040a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 16319b94acceSBarry Smith 16329665c990SLois Curfman McInnes Options Database Keys: 1633c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1634c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1635c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1636c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1637c7afd0dbSLois Curfman McInnes been hardwired into a code by 1638c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1639c7afd0dbSLois Curfman McInnes does not cancel those set via 1640c7afd0dbSLois Curfman McInnes the options database. 16419665c990SLois Curfman McInnes 1642639f9d9dSBarry Smith Notes: 16436bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 16446bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 16456bc08f3fSLois Curfman McInnes order in which they were set. 1646639f9d9dSBarry Smith 164736851e7fSLois Curfman McInnes Level: intermediate 164836851e7fSLois Curfman McInnes 16499b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 16509b94acceSBarry Smith 16515cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 16529b94acceSBarry Smith @*/ 165374679c65SBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx ) 16549b94acceSBarry Smith { 16553a40ed3dSBarry Smith PetscFunctionBegin; 1656639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1657a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set"); 1658639f9d9dSBarry Smith } 1659639f9d9dSBarry Smith 1660639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1661639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 16623a40ed3dSBarry Smith PetscFunctionReturn(0); 16639b94acceSBarry Smith } 16649b94acceSBarry Smith 16655615d1e5SSatish Balay #undef __FUNC__ 16665cd90555SBarry Smith #define __FUNC__ "SNESClearMonitor" 16675cd90555SBarry Smith /*@C 16685cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 16695cd90555SBarry Smith 1670c7afd0dbSLois Curfman McInnes Collective on SNES 1671c7afd0dbSLois Curfman McInnes 16725cd90555SBarry Smith Input Parameters: 16735cd90555SBarry Smith . snes - the SNES context 16745cd90555SBarry Smith 16755cd90555SBarry Smith Options Database: 1676c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1677c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1678c7afd0dbSLois Curfman McInnes set via the options database 16795cd90555SBarry Smith 16805cd90555SBarry Smith Notes: 16815cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 16825cd90555SBarry Smith 168336851e7fSLois Curfman McInnes Level: intermediate 168436851e7fSLois Curfman McInnes 16855cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 16865cd90555SBarry Smith 16875cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 16885cd90555SBarry Smith @*/ 16895cd90555SBarry Smith int SNESClearMonitor( SNES snes ) 16905cd90555SBarry Smith { 16915cd90555SBarry Smith PetscFunctionBegin; 16925cd90555SBarry Smith snes->numbermonitors = 0; 16935cd90555SBarry Smith PetscFunctionReturn(0); 16945cd90555SBarry Smith } 16955cd90555SBarry Smith 16965cd90555SBarry Smith #undef __FUNC__ 1697d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest" 16989b94acceSBarry Smith /*@C 16999b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 17009b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 17019b94acceSBarry Smith 1702fee21e36SBarry Smith Collective on SNES 1703fee21e36SBarry Smith 1704c7afd0dbSLois Curfman McInnes Input Parameters: 1705c7afd0dbSLois Curfman McInnes + snes - the SNES context 1706c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1707c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1708c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 17099b94acceSBarry Smith 1710c7afd0dbSLois Curfman McInnes Calling sequence of func: 17118d76a1e5SLois Curfman McInnes $ int func (SNES snes,double xnorm,double gnorm,double f,void *cctx) 1712c7afd0dbSLois Curfman McInnes 1713c7afd0dbSLois Curfman McInnes + snes - the SNES context 1714c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1715c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 1716c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods) 1717c7afd0dbSLois Curfman McInnes . f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods) 1718c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods) 1719c7afd0dbSLois Curfman McInnes - f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods) 17209b94acceSBarry Smith 172136851e7fSLois Curfman McInnes Level: advanced 172236851e7fSLois Curfman McInnes 17239b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 17249b94acceSBarry Smith 172540191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 172640191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 17279b94acceSBarry Smith @*/ 172874679c65SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,void*),void *cctx) 17299b94acceSBarry Smith { 17303a40ed3dSBarry Smith PetscFunctionBegin; 17319b94acceSBarry Smith (snes)->converged = func; 17329b94acceSBarry Smith (snes)->cnvP = cctx; 17333a40ed3dSBarry Smith PetscFunctionReturn(0); 17349b94acceSBarry Smith } 17359b94acceSBarry Smith 17365615d1e5SSatish Balay #undef __FUNC__ 1737d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" 1738c9005455SLois Curfman McInnes /*@ 1739c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1740c9005455SLois Curfman McInnes 1741fee21e36SBarry Smith Collective on SNES 1742fee21e36SBarry Smith 1743c7afd0dbSLois Curfman McInnes Input Parameters: 1744c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1745c7afd0dbSLois Curfman McInnes . a - array to hold history 1746758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1747758f92a0SBarry Smith negative if not converged) for each solve. 1748758f92a0SBarry Smith . na - size of a and its 1749758f92a0SBarry Smith - reset - PETSC_TRUTH indicates each new nonlinear solve resets the history counter to zero, 1750758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1751c7afd0dbSLois Curfman McInnes 1752c9005455SLois Curfman McInnes Notes: 1753c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1754c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1755c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1756c9005455SLois Curfman McInnes at each step. 1757c9005455SLois Curfman McInnes 1758c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1759c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1760c9005455SLois Curfman McInnes during the section of code that is being timed. 1761c9005455SLois Curfman McInnes 176236851e7fSLois Curfman McInnes Level: intermediate 176336851e7fSLois Curfman McInnes 1764c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1765758f92a0SBarry Smith 176608405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1767758f92a0SBarry Smith 1768c9005455SLois Curfman McInnes @*/ 1769758f92a0SBarry Smith int SNESSetConvergenceHistory(SNES snes, double *a, int *its,int na,PetscTruth reset) 1770c9005455SLois Curfman McInnes { 17713a40ed3dSBarry Smith PetscFunctionBegin; 1772c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1773c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1774c9005455SLois Curfman McInnes snes->conv_hist = a; 1775758f92a0SBarry Smith snes->conv_hist_its = its; 1776758f92a0SBarry Smith snes->conv_hist_max = na; 1777758f92a0SBarry Smith snes->conv_hist_reset = reset; 1778758f92a0SBarry Smith PetscFunctionReturn(0); 1779758f92a0SBarry Smith } 1780758f92a0SBarry Smith 1781758f92a0SBarry Smith #undef __FUNC__ 1782758f92a0SBarry Smith #define __FUNC__ "SNESGetConvergenceHistory" 17830c4c9dddSBarry Smith /*@C 1784758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1785758f92a0SBarry Smith 1786758f92a0SBarry Smith Collective on SNES 1787758f92a0SBarry Smith 1788758f92a0SBarry Smith Input Parameter: 1789758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1790758f92a0SBarry Smith 1791758f92a0SBarry Smith Output Parameters: 1792758f92a0SBarry Smith . a - array to hold history 1793758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1794758f92a0SBarry Smith negative if not converged) for each solve. 1795758f92a0SBarry Smith - na - size of a and its 1796758f92a0SBarry Smith 1797758f92a0SBarry Smith Notes: 1798758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1799758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1800758f92a0SBarry Smith 1801758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1802758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1803758f92a0SBarry Smith during the section of code that is being timed. 1804758f92a0SBarry Smith 1805758f92a0SBarry Smith Level: intermediate 1806758f92a0SBarry Smith 1807758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1808758f92a0SBarry Smith 1809758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1810758f92a0SBarry Smith 1811758f92a0SBarry Smith @*/ 1812758f92a0SBarry Smith int SNESGetConvergenceHistory(SNES snes, double **a, int **its,int *na) 1813758f92a0SBarry Smith { 1814758f92a0SBarry Smith PetscFunctionBegin; 1815758f92a0SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1816758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1817758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1818758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 18193a40ed3dSBarry Smith PetscFunctionReturn(0); 1820c9005455SLois Curfman McInnes } 1821c9005455SLois Curfman McInnes 1822c9005455SLois Curfman McInnes #undef __FUNC__ 18235615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 18249b94acceSBarry Smith /* 18259b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 18269b94acceSBarry Smith positive parameter delta. 18279b94acceSBarry Smith 18289b94acceSBarry Smith Input Parameters: 1829c7afd0dbSLois Curfman McInnes + snes - the SNES context 18309b94acceSBarry Smith . y - approximate solution of linear system 18319b94acceSBarry Smith . fnorm - 2-norm of current function 1832c7afd0dbSLois Curfman McInnes - delta - trust region size 18339b94acceSBarry Smith 18349b94acceSBarry Smith Output Parameters: 1835c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 18369b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 18379b94acceSBarry Smith region, and exceeds zero otherwise. 1838c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 18399b94acceSBarry Smith 18409b94acceSBarry Smith Note: 184140191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 18429b94acceSBarry Smith is set to be the maximum allowable step size. 18439b94acceSBarry Smith 18449b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 18459b94acceSBarry Smith */ 18469b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 18479b94acceSBarry Smith double *gpnorm,double *ynorm) 18489b94acceSBarry Smith { 18499b94acceSBarry Smith double norm; 18509b94acceSBarry Smith Scalar cnorm; 18513a40ed3dSBarry Smith int ierr; 18523a40ed3dSBarry Smith 18533a40ed3dSBarry Smith PetscFunctionBegin; 18543a40ed3dSBarry Smith ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr); 18559b94acceSBarry Smith if (norm > *delta) { 18569b94acceSBarry Smith norm = *delta/norm; 18579b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 18589b94acceSBarry Smith cnorm = norm; 18599b94acceSBarry Smith VecScale( &cnorm, y ); 18609b94acceSBarry Smith *ynorm = *delta; 18619b94acceSBarry Smith } else { 18629b94acceSBarry Smith *gpnorm = 0.0; 18639b94acceSBarry Smith *ynorm = norm; 18649b94acceSBarry Smith } 18653a40ed3dSBarry Smith PetscFunctionReturn(0); 18669b94acceSBarry Smith } 18679b94acceSBarry Smith 18685615d1e5SSatish Balay #undef __FUNC__ 18695615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 18709b94acceSBarry Smith /*@ 18719b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 187263a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 18739b94acceSBarry Smith 1874c7afd0dbSLois Curfman McInnes Collective on SNES 1875c7afd0dbSLois Curfman McInnes 1876b2002411SLois Curfman McInnes Input Parameters: 1877c7afd0dbSLois Curfman McInnes + snes - the SNES context 1878c7afd0dbSLois Curfman McInnes - x - the solution vector 18799b94acceSBarry Smith 18809b94acceSBarry Smith Output Parameter: 1881b2002411SLois Curfman McInnes . its - number of iterations until termination 18829b94acceSBarry Smith 1883b2002411SLois Curfman McInnes Notes: 18848ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 18858ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 18868ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 18878ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 18888ddd3da0SLois Curfman McInnes 188936851e7fSLois Curfman McInnes Level: beginner 189036851e7fSLois Curfman McInnes 18919b94acceSBarry Smith .keywords: SNES, nonlinear, solve 18929b94acceSBarry Smith 189363a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 18949b94acceSBarry Smith @*/ 18958ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 18969b94acceSBarry Smith { 18973c7409f5SSatish Balay int ierr, flg; 1898052efed2SBarry Smith 18993a40ed3dSBarry Smith PetscFunctionBegin; 190077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 190174679c65SBarry Smith PetscValidIntPointer(its); 190282bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x); CHKERRQ(ierr);} 1903c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 1904758f92a0SBarry Smith if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 19059b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 1906c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 19079b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its); CHKERRQ(ierr); 19089b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 19093c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg); CHKERRQ(ierr); 19106d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD); CHKERRQ(ierr); } 19113a40ed3dSBarry Smith PetscFunctionReturn(0); 19129b94acceSBarry Smith } 19139b94acceSBarry Smith 19149b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 19159b94acceSBarry Smith 19165615d1e5SSatish Balay #undef __FUNC__ 19175615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 191882bf6240SBarry Smith /*@C 19194b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 19209b94acceSBarry Smith 1921fee21e36SBarry Smith Collective on SNES 1922fee21e36SBarry Smith 1923c7afd0dbSLois Curfman McInnes Input Parameters: 1924c7afd0dbSLois Curfman McInnes + snes - the SNES context 1925c7afd0dbSLois Curfman McInnes - method - a known method 1926c7afd0dbSLois Curfman McInnes 1927c7afd0dbSLois Curfman McInnes Options Database Key: 1928c7afd0dbSLois Curfman McInnes . -snes_type <method> - Sets the method; use -help for a list 1929c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 1930ae12b187SLois Curfman McInnes 19319b94acceSBarry Smith Notes: 19329b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 193336851e7fSLois Curfman McInnes + SNES_EQ_LS - Newton's method with line search 1934c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 1935c7afd0dbSLois Curfman McInnes . SNES_EQ_TR - Newton's method with trust region 1936c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 1937c7afd0dbSLois Curfman McInnes . SNES_UM_TR - Newton's method with trust region 1938c7afd0dbSLois Curfman McInnes (unconstrained minimization) 193936851e7fSLois Curfman McInnes - SNES_UM_LS - Newton's method with line search 1940c7afd0dbSLois Curfman McInnes (unconstrained minimization) 19419b94acceSBarry Smith 1942ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 1943ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 1944ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 1945ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 1946ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 1947ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 1948ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 1949ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 1950ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 195136851e7fSLois Curfman McInnes appropriate method. In other words, this routine is not for beginners. 195236851e7fSLois Curfman McInnes 195336851e7fSLois Curfman McInnes Level: intermediate 1954a703fe33SLois Curfman McInnes 1955f536c699SSatish Balay .keywords: SNES, set, method 19569b94acceSBarry Smith @*/ 19574b0e389bSBarry Smith int SNESSetType(SNES snes,SNESType method) 19589b94acceSBarry Smith { 195984cb2905SBarry Smith int ierr; 19609b94acceSBarry Smith int (*r)(SNES); 19613a40ed3dSBarry Smith 19623a40ed3dSBarry Smith PetscFunctionBegin; 196377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 196482bf6240SBarry Smith 19653f1db9ecSBarry Smith if (PetscTypeCompare(snes->type_name,method)) PetscFunctionReturn(0); 196682bf6240SBarry Smith 196782bf6240SBarry Smith if (snes->setupcalled) { 1968e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 196982bf6240SBarry Smith snes->data = 0; 197082bf6240SBarry Smith } 19717d1a2b2bSBarry Smith 19729b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 197382bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);} 197482bf6240SBarry Smith 1975488ecbafSBarry Smith ierr = FListFind(snes->comm, SNESList, method,(int (**)(void *)) &r );CHKERRQ(ierr); 197682bf6240SBarry Smith 1977596552b5SBarry Smith if (!r) SETERRQ1(1,1,"Unable to find requested SNES type %s",method); 19781548b14aSBarry Smith 19790452661fSBarry Smith if (snes->data) PetscFree(snes->data); 198082bf6240SBarry Smith snes->data = 0; 19813a40ed3dSBarry Smith ierr = (*r)(snes); CHKERRQ(ierr); 198282bf6240SBarry Smith 198382bf6240SBarry Smith if (snes->type_name) PetscFree(snes->type_name); 198482bf6240SBarry Smith snes->type_name = (char *) PetscMalloc((PetscStrlen(method)+1)*sizeof(char));CHKPTRQ(snes->type_name); 198582bf6240SBarry Smith PetscStrcpy(snes->type_name,method); 198682bf6240SBarry Smith snes->set_method_called = 1; 198782bf6240SBarry Smith 19883a40ed3dSBarry Smith PetscFunctionReturn(0); 19899b94acceSBarry Smith } 19909b94acceSBarry Smith 1991a847f771SSatish Balay 19929b94acceSBarry Smith /* --------------------------------------------------------------------- */ 19935615d1e5SSatish Balay #undef __FUNC__ 1994d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy" 19959b94acceSBarry Smith /*@C 19969b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 19979b94acceSBarry Smith registered by SNESRegister(). 19989b94acceSBarry Smith 1999fee21e36SBarry Smith Not Collective 2000fee21e36SBarry Smith 200136851e7fSLois Curfman McInnes Level: advanced 200236851e7fSLois Curfman McInnes 20039b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 20049b94acceSBarry Smith 20059b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 20069b94acceSBarry Smith @*/ 2007cf256101SBarry Smith int SNESRegisterDestroy(void) 20089b94acceSBarry Smith { 200982bf6240SBarry Smith int ierr; 201082bf6240SBarry Smith 20113a40ed3dSBarry Smith PetscFunctionBegin; 201282bf6240SBarry Smith if (SNESList) { 2013488ecbafSBarry Smith ierr = FListDestroy( SNESList );CHKERRQ(ierr); 201482bf6240SBarry Smith SNESList = 0; 20159b94acceSBarry Smith } 201684cb2905SBarry Smith SNESRegisterAllCalled = 0; 20173a40ed3dSBarry Smith PetscFunctionReturn(0); 20189b94acceSBarry Smith } 20199b94acceSBarry Smith 20205615d1e5SSatish Balay #undef __FUNC__ 2021d4bb536fSBarry Smith #define __FUNC__ "SNESGetType" 20229b94acceSBarry Smith /*@C 20239a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 20249b94acceSBarry Smith 2025c7afd0dbSLois Curfman McInnes Not Collective 2026c7afd0dbSLois Curfman McInnes 20279b94acceSBarry Smith Input Parameter: 20284b0e389bSBarry Smith . snes - nonlinear solver context 20299b94acceSBarry Smith 20309b94acceSBarry Smith Output Parameter: 203182bf6240SBarry Smith . method - SNES method (a charactor string) 20329b94acceSBarry Smith 203336851e7fSLois Curfman McInnes Level: intermediate 203436851e7fSLois Curfman McInnes 20359b94acceSBarry Smith .keywords: SNES, nonlinear, get, method, name 20369b94acceSBarry Smith @*/ 203782bf6240SBarry Smith int SNESGetType(SNES snes, SNESType *method) 20389b94acceSBarry Smith { 20393a40ed3dSBarry Smith PetscFunctionBegin; 204082bf6240SBarry Smith *method = snes->type_name; 20413a40ed3dSBarry Smith PetscFunctionReturn(0); 20429b94acceSBarry Smith } 20439b94acceSBarry Smith 20445615d1e5SSatish Balay #undef __FUNC__ 2045d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution" 20469b94acceSBarry Smith /*@C 20479b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 20489b94acceSBarry Smith stored. 20499b94acceSBarry Smith 2050c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2051c7afd0dbSLois Curfman McInnes 20529b94acceSBarry Smith Input Parameter: 20539b94acceSBarry Smith . snes - the SNES context 20549b94acceSBarry Smith 20559b94acceSBarry Smith Output Parameter: 20569b94acceSBarry Smith . x - the solution 20579b94acceSBarry Smith 205836851e7fSLois Curfman McInnes Level: advanced 205936851e7fSLois Curfman McInnes 20609b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 20619b94acceSBarry Smith 2062a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 20639b94acceSBarry Smith @*/ 20649b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 20659b94acceSBarry Smith { 20663a40ed3dSBarry Smith PetscFunctionBegin; 206777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 20689b94acceSBarry Smith *x = snes->vec_sol_always; 20693a40ed3dSBarry Smith PetscFunctionReturn(0); 20709b94acceSBarry Smith } 20719b94acceSBarry Smith 20725615d1e5SSatish Balay #undef __FUNC__ 2073d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" 20749b94acceSBarry Smith /*@C 20759b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 20769b94acceSBarry Smith stored. 20779b94acceSBarry Smith 2078c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2079c7afd0dbSLois Curfman McInnes 20809b94acceSBarry Smith Input Parameter: 20819b94acceSBarry Smith . snes - the SNES context 20829b94acceSBarry Smith 20839b94acceSBarry Smith Output Parameter: 20849b94acceSBarry Smith . x - the solution update 20859b94acceSBarry Smith 208636851e7fSLois Curfman McInnes Level: advanced 208736851e7fSLois Curfman McInnes 20889b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 20899b94acceSBarry Smith 20909b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 20919b94acceSBarry Smith @*/ 20929b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 20939b94acceSBarry Smith { 20943a40ed3dSBarry Smith PetscFunctionBegin; 209577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 20969b94acceSBarry Smith *x = snes->vec_sol_update_always; 20973a40ed3dSBarry Smith PetscFunctionReturn(0); 20989b94acceSBarry Smith } 20999b94acceSBarry Smith 21005615d1e5SSatish Balay #undef __FUNC__ 2101d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction" 21029b94acceSBarry Smith /*@C 21033638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 21049b94acceSBarry Smith 2105c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2106c7afd0dbSLois Curfman McInnes 21079b94acceSBarry Smith Input Parameter: 21089b94acceSBarry Smith . snes - the SNES context 21099b94acceSBarry Smith 21109b94acceSBarry Smith Output Parameter: 21113638b69dSLois Curfman McInnes . r - the function 21129b94acceSBarry Smith 21139b94acceSBarry Smith Notes: 21149b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 21159b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 21169b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 21179b94acceSBarry Smith 211836851e7fSLois Curfman McInnes Level: advanced 211936851e7fSLois Curfman McInnes 2120a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 21219b94acceSBarry Smith 212231615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 212331615d04SLois Curfman McInnes SNESGetGradient() 21249b94acceSBarry Smith @*/ 21259b94acceSBarry Smith int SNESGetFunction(SNES snes,Vec *r) 21269b94acceSBarry Smith { 21273a40ed3dSBarry Smith PetscFunctionBegin; 212877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2129a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 2130a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 2131a8c6a408SBarry Smith } 21329b94acceSBarry Smith *r = snes->vec_func_always; 21333a40ed3dSBarry Smith PetscFunctionReturn(0); 21349b94acceSBarry Smith } 21359b94acceSBarry Smith 21365615d1e5SSatish Balay #undef __FUNC__ 2137d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient" 21389b94acceSBarry Smith /*@C 21393638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 21409b94acceSBarry Smith 2141c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2142c7afd0dbSLois Curfman McInnes 21439b94acceSBarry Smith Input Parameter: 21449b94acceSBarry Smith . snes - the SNES context 21459b94acceSBarry Smith 21469b94acceSBarry Smith Output Parameter: 21479b94acceSBarry Smith . r - the gradient 21489b94acceSBarry Smith 21499b94acceSBarry Smith Notes: 21509b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 21519b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 21529b94acceSBarry Smith SNESGetFunction(). 21539b94acceSBarry Smith 215436851e7fSLois Curfman McInnes Level: advanced 215536851e7fSLois Curfman McInnes 21569b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 21579b94acceSBarry Smith 215831615d04SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction() 21599b94acceSBarry Smith @*/ 21609b94acceSBarry Smith int SNESGetGradient(SNES snes,Vec *r) 21619b94acceSBarry Smith { 21623a40ed3dSBarry Smith PetscFunctionBegin; 216377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 21643a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 2165a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 21663a40ed3dSBarry Smith } 21679b94acceSBarry Smith *r = snes->vec_func_always; 21683a40ed3dSBarry Smith PetscFunctionReturn(0); 21699b94acceSBarry Smith } 21709b94acceSBarry Smith 21715615d1e5SSatish Balay #undef __FUNC__ 2172d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" 21739b94acceSBarry Smith /*@ 21749b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 21759b94acceSBarry Smith unconstrained minimization problems. 21769b94acceSBarry Smith 2177c7afd0dbSLois Curfman McInnes Not Collective 2178c7afd0dbSLois Curfman McInnes 21799b94acceSBarry Smith Input Parameter: 21809b94acceSBarry Smith . snes - the SNES context 21819b94acceSBarry Smith 21829b94acceSBarry Smith Output Parameter: 21839b94acceSBarry Smith . r - the function 21849b94acceSBarry Smith 21859b94acceSBarry Smith Notes: 21869b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 21879b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 21889b94acceSBarry Smith SNESGetFunction(). 21899b94acceSBarry Smith 219036851e7fSLois Curfman McInnes Level: advanced 219136851e7fSLois Curfman McInnes 21929b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 21939b94acceSBarry Smith 219431615d04SLois Curfman McInnes .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction() 21959b94acceSBarry Smith @*/ 21969b94acceSBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r) 21979b94acceSBarry Smith { 21983a40ed3dSBarry Smith PetscFunctionBegin; 219977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 220074679c65SBarry Smith PetscValidScalarPointer(r); 22013a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 2202a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 22033a40ed3dSBarry Smith } 22049b94acceSBarry Smith *r = snes->fc; 22053a40ed3dSBarry Smith PetscFunctionReturn(0); 22069b94acceSBarry Smith } 22079b94acceSBarry Smith 22085615d1e5SSatish Balay #undef __FUNC__ 2209d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" 22103c7409f5SSatish Balay /*@C 22113c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 2212d850072dSLois Curfman McInnes SNES options in the database. 22133c7409f5SSatish Balay 2214fee21e36SBarry Smith Collective on SNES 2215fee21e36SBarry Smith 2216c7afd0dbSLois Curfman McInnes Input Parameter: 2217c7afd0dbSLois Curfman McInnes + snes - the SNES context 2218c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2219c7afd0dbSLois Curfman McInnes 2220d850072dSLois Curfman McInnes Notes: 2221a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2222c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2223d850072dSLois Curfman McInnes 222436851e7fSLois Curfman McInnes Level: advanced 222536851e7fSLois Curfman McInnes 22263c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 2227a86d99e1SLois Curfman McInnes 2228a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 22293c7409f5SSatish Balay @*/ 22303c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 22313c7409f5SSatish Balay { 22323c7409f5SSatish Balay int ierr; 22333c7409f5SSatish Balay 22343a40ed3dSBarry Smith PetscFunctionBegin; 223577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2236639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 22373c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 22383a40ed3dSBarry Smith PetscFunctionReturn(0); 22393c7409f5SSatish Balay } 22403c7409f5SSatish Balay 22415615d1e5SSatish Balay #undef __FUNC__ 2242d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" 22433c7409f5SSatish Balay /*@C 2244f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2245d850072dSLois Curfman McInnes SNES options in the database. 22463c7409f5SSatish Balay 2247fee21e36SBarry Smith Collective on SNES 2248fee21e36SBarry Smith 2249c7afd0dbSLois Curfman McInnes Input Parameters: 2250c7afd0dbSLois Curfman McInnes + snes - the SNES context 2251c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2252c7afd0dbSLois Curfman McInnes 2253d850072dSLois Curfman McInnes Notes: 2254a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2255c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2256d850072dSLois Curfman McInnes 225736851e7fSLois Curfman McInnes Level: advanced 225836851e7fSLois Curfman McInnes 22593c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2260a86d99e1SLois Curfman McInnes 2261a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 22623c7409f5SSatish Balay @*/ 22633c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 22643c7409f5SSatish Balay { 22653c7409f5SSatish Balay int ierr; 22663c7409f5SSatish Balay 22673a40ed3dSBarry Smith PetscFunctionBegin; 226877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2269639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 22703c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 22713a40ed3dSBarry Smith PetscFunctionReturn(0); 22723c7409f5SSatish Balay } 22733c7409f5SSatish Balay 22745615d1e5SSatish Balay #undef __FUNC__ 2275d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" 2276*9ab63eb5SSatish Balay /*@C 22773c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 22783c7409f5SSatish Balay SNES options in the database. 22793c7409f5SSatish Balay 2280c7afd0dbSLois Curfman McInnes Not Collective 2281c7afd0dbSLois Curfman McInnes 22823c7409f5SSatish Balay Input Parameter: 22833c7409f5SSatish Balay . snes - the SNES context 22843c7409f5SSatish Balay 22853c7409f5SSatish Balay Output Parameter: 22863c7409f5SSatish Balay . prefix - pointer to the prefix string used 22873c7409f5SSatish Balay 2288*9ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 2289*9ab63eb5SSatish Balay sufficient length to hold the prefix. 2290*9ab63eb5SSatish Balay 229136851e7fSLois Curfman McInnes Level: advanced 229236851e7fSLois Curfman McInnes 22933c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2294a86d99e1SLois Curfman McInnes 2295a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 22963c7409f5SSatish Balay @*/ 22973c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 22983c7409f5SSatish Balay { 22993c7409f5SSatish Balay int ierr; 23003c7409f5SSatish Balay 23013a40ed3dSBarry Smith PetscFunctionBegin; 230277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2303639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix); CHKERRQ(ierr); 23043a40ed3dSBarry Smith PetscFunctionReturn(0); 23053c7409f5SSatish Balay } 23063c7409f5SSatish Balay 2307ca161407SBarry Smith #undef __FUNC__ 2308ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp" 2309ca161407SBarry Smith /*@ 2310ca161407SBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2311ca161407SBarry Smith 2312c7afd0dbSLois Curfman McInnes Collective on SNES 2313c7afd0dbSLois Curfman McInnes 2314ca161407SBarry Smith Input Parameter: 2315ca161407SBarry Smith . snes - the SNES context 2316ca161407SBarry Smith 2317ca161407SBarry Smith Options Database Keys: 2318c7afd0dbSLois Curfman McInnes + -help - Prints SNES options 2319c7afd0dbSLois Curfman McInnes - -h - Prints SNES options 2320ca161407SBarry Smith 232136851e7fSLois Curfman McInnes Level: beginner 232236851e7fSLois Curfman McInnes 2323ca161407SBarry Smith .keywords: SNES, nonlinear, help 2324ca161407SBarry Smith 2325ca161407SBarry Smith .seealso: SNESSetFromOptions() 2326ca161407SBarry Smith @*/ 2327ca161407SBarry Smith int SNESPrintHelp(SNES snes) 2328ca161407SBarry Smith { 2329ca161407SBarry Smith char p[64]; 2330ca161407SBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2331ca161407SBarry Smith int ierr; 2332ca161407SBarry Smith 2333ca161407SBarry Smith PetscFunctionBegin; 2334ca161407SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2335ca161407SBarry Smith 2336ca161407SBarry Smith PetscStrcpy(p,"-"); 2337ca161407SBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2338ca161407SBarry Smith 2339ca161407SBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2340ca161407SBarry Smith 2341ebb8b11fSBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL); CHKERRQ(ierr);} 234276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n"); 2343488ecbafSBarry Smith ierr = FListPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr); 234476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p); 234576be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its); 234676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs); 234776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol); 234876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol); 234976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol); 235076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol); 235176be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n"); 235276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p); 235376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 2354ca161407SBarry Smith residual norm at each iteration.\n",p); 235576be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 2356ca161407SBarry Smith residual norm for small residual norms. This is useful to conceal\n\ 2357ca161407SBarry Smith meaningless digits that may be different on different machines.\n",p); 235876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p); 2359ca161407SBarry Smith if (snes->type == SNES_NONLINEAR_EQUATIONS) { 236076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2361ca161407SBarry Smith " Options for solving systems of nonlinear equations only:\n"); 236276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p); 236376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p); 236476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p); 2365e24b481bSBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf_ksp_monitor - if using matrix-free multiply then prints h at each KSP iteration\n",p); 236676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p); 236776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p); 236876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2369ca161407SBarry Smith " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version); 237076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2371ca161407SBarry Smith " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0); 237276be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2373ca161407SBarry Smith " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max); 237476be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2375ca161407SBarry Smith " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma); 237676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2377ca161407SBarry Smith " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha); 237876be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2379ca161407SBarry Smith " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2); 238076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2381ca161407SBarry Smith " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold); 2382ca161407SBarry Smith } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 238376be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm, 2384ca161407SBarry Smith " Options for solving unconstrained minimization problems only:\n"); 238576be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin); 238676be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p); 238776be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p); 2388ca161407SBarry Smith } 238976be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <method> for help on ",p); 239076be9ce4SBarry Smith (*PetscHelpPrintf)(snes->comm,"a particular method\n"); 2391ca161407SBarry Smith if (snes->printhelp) { 2392ca161407SBarry Smith ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr); 2393ca161407SBarry Smith } 2394ca161407SBarry Smith PetscFunctionReturn(0); 2395ca161407SBarry Smith } 23963c7409f5SSatish Balay 2397acb85ae6SSatish Balay /*MC 2398b2002411SLois Curfman McInnes SNESRegister - Adds a method to the nonlinear solver package. 23999b94acceSBarry Smith 2400b2002411SLois Curfman McInnes Synopsis: 2401b2002411SLois Curfman McInnes SNESRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES)) 24029b94acceSBarry Smith 24038d76a1e5SLois Curfman McInnes Not collective 24048d76a1e5SLois Curfman McInnes 2405b2002411SLois Curfman McInnes Input Parameters: 2406c7afd0dbSLois Curfman McInnes + name_solver - name of a new user-defined solver 2407b2002411SLois Curfman McInnes . path - path (either absolute or relative) the library containing this solver 2408b2002411SLois Curfman McInnes . name_create - name of routine to create method context 2409c7afd0dbSLois Curfman McInnes - routine_create - routine to create method context 24109b94acceSBarry Smith 2411b2002411SLois Curfman McInnes Notes: 2412b2002411SLois Curfman McInnes SNESRegister() may be called multiple times to add several user-defined solvers. 24133a40ed3dSBarry Smith 2414b2002411SLois Curfman McInnes If dynamic libraries are used, then the fourth input argument (routine_create) 2415b2002411SLois Curfman McInnes is ignored. 2416b2002411SLois Curfman McInnes 2417b2002411SLois Curfman McInnes Sample usage: 241869225d78SLois Curfman McInnes .vb 2419b2002411SLois Curfman McInnes SNESRegister("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a, 2420b2002411SLois Curfman McInnes "MySolverCreate",MySolverCreate); 242169225d78SLois Curfman McInnes .ve 2422b2002411SLois Curfman McInnes 2423b2002411SLois Curfman McInnes Then, your solver can be chosen with the procedural interface via 2424b2002411SLois Curfman McInnes $ SNESSetType(snes,"my_solver") 2425b2002411SLois Curfman McInnes or at runtime via the option 2426b2002411SLois Curfman McInnes $ -snes_type my_solver 2427b2002411SLois Curfman McInnes 242836851e7fSLois Curfman McInnes Level: advanced 242936851e7fSLois Curfman McInnes 2430dd438238SBarry Smith $PETSC_ARCH and $BOPT occuring in pathname will be replaced with appropriate values. 2431dd438238SBarry Smith 2432b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register 2433b2002411SLois Curfman McInnes 2434b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy() 2435b2002411SLois Curfman McInnes M*/ 2436b2002411SLois Curfman McInnes 2437b2002411SLois Curfman McInnes #undef __FUNC__ 2438b2002411SLois Curfman McInnes #define __FUNC__ "SNESRegister_Private" 2439b2002411SLois Curfman McInnes int SNESRegister_Private(char *sname,char *path,char *name,int (*function)(SNES)) 2440b2002411SLois Curfman McInnes { 2441b2002411SLois Curfman McInnes char fullname[256]; 2442b2002411SLois Curfman McInnes int ierr; 2443b2002411SLois Curfman McInnes 2444b2002411SLois Curfman McInnes PetscFunctionBegin; 2445b2002411SLois Curfman McInnes PetscStrcpy(fullname,path); PetscStrcat(fullname,":");PetscStrcat(fullname,name); 2446488ecbafSBarry Smith ierr = FListAdd_Private(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr); 2447b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2448b2002411SLois Curfman McInnes } 2449