1a5eb4965SSatish Balay #ifdef PETSC_RCS_HEADER 2*454a90a3SBarry Smith static char vcid[] = "$Id: snes.c,v 1.196 1999/09/27 21:31:38 bsmith Exp bsmith $"; 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; 46*454a90a3SBarry Smith char *type; 479b94acceSBarry Smith 483a40ed3dSBarry Smith PetscFunctionBegin; 4974679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5074679c65SBarry Smith if (viewer) {PetscValidHeader(viewer);} 5174679c65SBarry Smith else { viewer = VIEWER_STDOUT_SELF; } 5274679c65SBarry Smith 53*454a90a3SBarry Smith if (PetscTypeCompare(viewer,ASCII_VIEWER)) { 54d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 55*454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 56*454a90a3SBarry Smith if (type) { 57*454a90a3SBarry Smith ierr = ViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 58184914b5SBarry Smith } else { 59*454a90a3SBarry Smith ierr = ViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 60184914b5SBarry Smith } 610ef38995SBarry Smith if (snes->view) { 620ef38995SBarry Smith ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 630ef38995SBarry Smith ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 640ef38995SBarry Smith ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 650ef38995SBarry Smith } 66d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 67d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 68d132466eSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol);CHKERRQ(ierr); 69d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr); 70d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr); 710ef38995SBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 72d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," min function tolerance=%g\n",snes->fmin);CHKERRQ(ierr); 730ef38995SBarry Smith } 749b94acceSBarry Smith if (snes->ksp_ewconv) { 759b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 769b94acceSBarry Smith if (kctx) { 77d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr); 78d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 79d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 809b94acceSBarry Smith } 819b94acceSBarry Smith } 82*454a90a3SBarry Smith } else if (PetscTypeCompare(viewer,STRING_VIEWER)) { 83*454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 84*454a90a3SBarry Smith ierr = ViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 8519bcc07fSBarry Smith } 8677ed5343SBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 870ef38995SBarry Smith ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 889b94acceSBarry Smith ierr = SLESView(sles,viewer);CHKERRQ(ierr); 890ef38995SBarry Smith ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 903a40ed3dSBarry Smith PetscFunctionReturn(0); 919b94acceSBarry Smith } 929b94acceSBarry Smith 93639f9d9dSBarry Smith /* 94639f9d9dSBarry Smith We retain a list of functions that also take SNES command 95639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 96639f9d9dSBarry Smith */ 97639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 98639f9d9dSBarry Smith static int numberofsetfromoptions; 99639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 100639f9d9dSBarry Smith 1015615d1e5SSatish Balay #undef __FUNC__ 102d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker" 103639f9d9dSBarry Smith /*@ 104639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 105639f9d9dSBarry Smith 106c7afd0dbSLois Curfman McInnes Not Collective 107c7afd0dbSLois Curfman McInnes 108639f9d9dSBarry Smith Input Parameter: 109639f9d9dSBarry Smith . snescheck - function that checks for options 110639f9d9dSBarry Smith 11136851e7fSLois Curfman McInnes Level: developer 11236851e7fSLois Curfman McInnes 113639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 114639f9d9dSBarry Smith @*/ 115639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 116639f9d9dSBarry Smith { 1173a40ed3dSBarry Smith PetscFunctionBegin; 118639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 119a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed"); 120639f9d9dSBarry Smith } 121639f9d9dSBarry Smith 122639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 1233a40ed3dSBarry Smith PetscFunctionReturn(0); 124639f9d9dSBarry Smith } 125639f9d9dSBarry Smith 1265615d1e5SSatish Balay #undef __FUNC__ 12715091d37SBarry Smith #define __FUNC__ "SNESSetTypeFromOptions" 12815091d37SBarry Smith /*@ 12915091d37SBarry Smith SNESSetTypeFromOptions - Sets the SNES solver type from the options database, 13015091d37SBarry Smith or sets a default if none is give. 13115091d37SBarry Smith 13215091d37SBarry Smith Collective on SNES 13315091d37SBarry Smith 13415091d37SBarry Smith Input Parameter: 13515091d37SBarry Smith . snes - the SNES context 13615091d37SBarry Smith 13715091d37SBarry Smith Options Database Keys: 13815091d37SBarry Smith . -snes_type <type> - SNES_EQ_LS, SNES_EQ_TR, SNES_UM_TR, SNES_UM_LS etc 13915091d37SBarry Smith 14015091d37SBarry Smith Level: beginner 14115091d37SBarry Smith 14215091d37SBarry Smith .keywords: SNES, nonlinear, set, options, database 14315091d37SBarry Smith 14415091d37SBarry Smith .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetFromOptions() 14515091d37SBarry Smith @*/ 14615091d37SBarry Smith int SNESSetTypeFromOptions(SNES snes) 14715091d37SBarry Smith { 148*454a90a3SBarry Smith char type[256]; 14915091d37SBarry Smith int ierr, flg; 15015091d37SBarry Smith 15115091d37SBarry Smith PetscFunctionBegin; 15215091d37SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 15315091d37SBarry Smith if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp()"); 154*454a90a3SBarry Smith ierr = OptionsGetString(snes->prefix,"-snes_type",type,256,&flg);CHKERRQ(ierr); 15515091d37SBarry Smith if (flg) { 156*454a90a3SBarry Smith ierr = SNESSetType(snes,(SNESType) type);CHKERRQ(ierr); 15715091d37SBarry Smith } 15815091d37SBarry Smith /* 15915091d37SBarry Smith If SNES type has not yet been set, set it now 16015091d37SBarry Smith */ 16115091d37SBarry Smith if (!snes->type_name) { 16215091d37SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 16315091d37SBarry Smith ierr = SNESSetType(snes,SNES_EQ_LS);CHKERRQ(ierr); 16415091d37SBarry Smith } else { 16515091d37SBarry Smith ierr = SNESSetType(snes,SNES_UM_TR);CHKERRQ(ierr); 16615091d37SBarry Smith } 16715091d37SBarry Smith } 16815091d37SBarry Smith PetscFunctionReturn(0); 16915091d37SBarry Smith } 17015091d37SBarry Smith 17115091d37SBarry Smith #undef __FUNC__ 1725615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1739b94acceSBarry Smith /*@ 174682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1759b94acceSBarry Smith 176c7afd0dbSLois Curfman McInnes Collective on SNES 177c7afd0dbSLois Curfman McInnes 1789b94acceSBarry Smith Input Parameter: 1799b94acceSBarry Smith . snes - the SNES context 1809b94acceSBarry Smith 18136851e7fSLois Curfman McInnes Options Database Keys: 182b39c3a46SLois Curfman McInnes + -snes_type <type> - SNES_EQ_LS, SNES_EQ_TR, SNES_UM_TR, SNES_UM_LS etc 18382738288SBarry Smith . -snes_stol - convergence tolerance in terms of the norm 18482738288SBarry Smith of the change in the solution between steps 185b39c3a46SLois Curfman McInnes . -snes_atol <atol> - absolute tolerance of residual norm 186b39c3a46SLois Curfman McInnes . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 187b39c3a46SLois Curfman McInnes . -snes_max_it <max_it> - maximum number of iterations 188b39c3a46SLois Curfman McInnes . -snes_max_funcs <max_funcs> - maximum number of function evaluations 189b39c3a46SLois Curfman McInnes . -snes_trtol <trtol> - trust region tolerance 19082738288SBarry Smith . -snes_no_convergence_test - skip convergence test in nonlinear or minimization 19182738288SBarry Smith solver; hence iterations will continue until max_it 1921fbbfb26SLois Curfman McInnes or some other criterion is reached. Saves expense 19382738288SBarry Smith of convergence test 19482738288SBarry Smith . -snes_monitor - prints residual norm at each iteration 195d132466eSBarry Smith . -snes_vecmonitor - plots solution at each iteration 196d132466eSBarry Smith . -snes_vecmonitor_update - plots update to solution at each iteration 19782738288SBarry Smith . -snes_xmonitor - plots residual norm at each iteration 198e24b481bSBarry Smith . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 19936851e7fSLois Curfman McInnes - -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 20082738288SBarry Smith 20182738288SBarry Smith Options Database for Eisenstat-Walker method: 20282738288SBarry Smith + -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence 20382738288SBarry Smith . -snes_ksp_eq_version ver - version of Eisenstat-Walker method 20436851e7fSLois Curfman McInnes . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 20536851e7fSLois Curfman McInnes . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 20636851e7fSLois Curfman McInnes . -snes_ksp_ew_gamma <gamma> - Sets gamma 20736851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha <alpha> - Sets alpha 20836851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 20936851e7fSLois Curfman McInnes - -snes_ksp_ew_threshold <threshold> - Sets threshold 21082738288SBarry Smith 21111ca99fdSLois Curfman McInnes Notes: 21211ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 21311ca99fdSLois Curfman McInnes the users manual. 21483e2fdc7SBarry Smith 21536851e7fSLois Curfman McInnes Level: beginner 21636851e7fSLois Curfman McInnes 2179b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 2189b94acceSBarry Smith 21915091d37SBarry Smith .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetTypeFromOptions() 2209b94acceSBarry Smith @*/ 2219b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 2229b94acceSBarry Smith { 2239b94acceSBarry Smith double tmp; 2249b94acceSBarry Smith SLES sles; 22581f57ec7SBarry Smith int ierr, flg,i,loc[4],nmax = 4; 2263c7409f5SSatish Balay int version = PETSC_DEFAULT; 2279b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 2289b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 2299b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 2309b94acceSBarry Smith double alpha = PETSC_DEFAULT; 2319b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 2329b94acceSBarry Smith double threshold = PETSC_DEFAULT; 2339b94acceSBarry Smith 2343a40ed3dSBarry Smith PetscFunctionBegin; 23577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 23615091d37SBarry Smith ierr = SNESSetTypeFromOptions(snes);CHKERRQ(ierr); 237ca161407SBarry Smith 23815091d37SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 2393c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg);CHKERRQ(ierr); 240d64ed03dSBarry Smith if (flg) { 241d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 242d64ed03dSBarry Smith } 2433c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg);CHKERRQ(ierr); 244d64ed03dSBarry Smith if (flg) { 245d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 246d64ed03dSBarry Smith } 2473c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);CHKERRQ(ierr); 248d64ed03dSBarry Smith if (flg) { 249d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 250d64ed03dSBarry Smith } 2513c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, &flg);CHKERRQ(ierr); 2523c7409f5SSatish Balay ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, &flg);CHKERRQ(ierr); 253d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg);CHKERRQ(ierr); 254d64ed03dSBarry Smith if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp);CHKERRQ(ierr); } 255d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg);CHKERRQ(ierr); 256d64ed03dSBarry Smith if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);CHKERRQ(ierr);} 2573c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg);CHKERRQ(ierr); 2583c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 259b18e04deSLois Curfman McInnes ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,&flg);CHKERRQ(ierr); 260b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,&flg);CHKERRQ(ierr); 261b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,&flg);CHKERRQ(ierr); 262b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,&flg);CHKERRQ(ierr); 263b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,&flg);CHKERRQ(ierr); 264b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,&flg);CHKERRQ(ierr); 265b18e04deSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,&flg);CHKERRQ(ierr); 26693c39befSBarry Smith 26793c39befSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg);CHKERRQ(ierr); 26893c39befSBarry Smith if (flg) {snes->converged = 0;} 26993c39befSBarry Smith 2709b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 2719b94acceSBarry Smith alpha2,threshold);CHKERRQ(ierr); 2725bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg);CHKERRQ(ierr); 2735cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 2743c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg);CHKERRQ(ierr); 275b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);} 2763c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg);CHKERRQ(ierr); 277b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);} 2783f1db9ecSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor",&flg);CHKERRQ(ierr); 279b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 280d132466eSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor_update",&flg);CHKERRQ(ierr); 281b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitorUpdate,0,0);CHKERRQ(ierr);} 28281f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2833c7409f5SSatish Balay if (flg) { 28417699dbbSLois Curfman McInnes int rank = 0; 285d7e8b826SBarry Smith DrawLG lg; 28617699dbbSLois Curfman McInnes MPI_Initialized(&rank); 287d132466eSBarry Smith if (rank) {ierr = MPI_Comm_rank(snes->comm,&rank);CHKERRQ(ierr);} 28817699dbbSLois Curfman McInnes if (!rank) { 28981f57ec7SBarry Smith ierr = SNESLGMonitorCreate(0,0,loc[0],loc[1],loc[2],loc[3],&lg);CHKERRQ(ierr); 290b8a78c4aSBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,lg,( int (*)(void *))SNESLGMonitorDestroy);CHKERRQ(ierr); 2919b94acceSBarry Smith PLogObjectParent(snes,lg); 2929b94acceSBarry Smith } 2939b94acceSBarry Smith } 294e24b481bSBarry Smith 2953c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_fd", &flg);CHKERRQ(ierr); 2963c7409f5SSatish Balay if (flg && snes->method_class == SNES_NONLINEAR_EQUATIONS) { 2979b94acceSBarry Smith ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre, 2989b94acceSBarry Smith SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 299981c4779SBarry Smith PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 300981c4779SBarry Smith } else if (flg && snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 30131615d04SLois Curfman McInnes ierr = SNESSetHessian(snes,snes->jacobian,snes->jacobian_pre, 30231615d04SLois Curfman McInnes SNESDefaultComputeHessian,snes->funP);CHKERRQ(ierr); 303d64ed03dSBarry Smith PLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Hessian matrix\n"); 3049b94acceSBarry Smith } 305639f9d9dSBarry Smith 306639f9d9dSBarry Smith for ( i=0; i<numberofsetfromoptions; i++ ) { 307639f9d9dSBarry Smith ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 308639f9d9dSBarry Smith } 309639f9d9dSBarry Smith 3109b94acceSBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 3119b94acceSBarry Smith ierr = SLESSetFromOptions(sles);CHKERRQ(ierr); 31293993e2dSLois Curfman McInnes 313e24b481bSBarry Smith /* set the special KSP monitor for matrix-free application */ 314e24b481bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_ksp_monitor",&flg);CHKERRQ(ierr); 315e24b481bSBarry Smith if (flg) { 316e24b481bSBarry Smith KSP ksp; 317e24b481bSBarry Smith ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 318b8a78c4aSBarry Smith ierr = KSPSetMonitor(ksp,MatSNESMFKSPMonitor,PETSC_NULL,0);CHKERRQ(ierr); 319e24b481bSBarry Smith } 320e24b481bSBarry Smith 32182bf6240SBarry Smith ierr = OptionsHasName(PETSC_NULL,"-help", &flg);CHKERRQ(ierr); 32282bf6240SBarry Smith if (flg) { ierr = SNESPrintHelp(snes);CHKERRQ(ierr);} 32382bf6240SBarry Smith 3243a40ed3dSBarry Smith if (snes->setfromoptions) { 3253a40ed3dSBarry Smith ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 3263a40ed3dSBarry Smith } 3273a40ed3dSBarry Smith PetscFunctionReturn(0); 3289b94acceSBarry Smith } 3299b94acceSBarry Smith 330a847f771SSatish Balay 3315615d1e5SSatish Balay #undef __FUNC__ 332d4bb536fSBarry Smith #define __FUNC__ "SNESSetApplicationContext" 3339b94acceSBarry Smith /*@ 3349b94acceSBarry Smith SNESSetApplicationContext - Sets the optional user-defined context for 3359b94acceSBarry Smith the nonlinear solvers. 3369b94acceSBarry Smith 337fee21e36SBarry Smith Collective on SNES 338fee21e36SBarry Smith 339c7afd0dbSLois Curfman McInnes Input Parameters: 340c7afd0dbSLois Curfman McInnes + snes - the SNES context 341c7afd0dbSLois Curfman McInnes - usrP - optional user context 342c7afd0dbSLois Curfman McInnes 34336851e7fSLois Curfman McInnes Level: intermediate 34436851e7fSLois Curfman McInnes 3459b94acceSBarry Smith .keywords: SNES, nonlinear, set, application, context 3469b94acceSBarry Smith 3479b94acceSBarry Smith .seealso: SNESGetApplicationContext() 3489b94acceSBarry Smith @*/ 3499b94acceSBarry Smith int SNESSetApplicationContext(SNES snes,void *usrP) 3509b94acceSBarry Smith { 3513a40ed3dSBarry Smith PetscFunctionBegin; 35277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3539b94acceSBarry Smith snes->user = usrP; 3543a40ed3dSBarry Smith PetscFunctionReturn(0); 3559b94acceSBarry Smith } 35674679c65SBarry Smith 3575615d1e5SSatish Balay #undef __FUNC__ 358d4bb536fSBarry Smith #define __FUNC__ "SNESGetApplicationContext" 3599b94acceSBarry Smith /*@C 3609b94acceSBarry Smith SNESGetApplicationContext - Gets the user-defined context for the 3619b94acceSBarry Smith nonlinear solvers. 3629b94acceSBarry Smith 363c7afd0dbSLois Curfman McInnes Not Collective 364c7afd0dbSLois Curfman McInnes 3659b94acceSBarry Smith Input Parameter: 3669b94acceSBarry Smith . snes - SNES context 3679b94acceSBarry Smith 3689b94acceSBarry Smith Output Parameter: 3699b94acceSBarry Smith . usrP - user context 3709b94acceSBarry Smith 37136851e7fSLois Curfman McInnes Level: intermediate 37236851e7fSLois Curfman McInnes 3739b94acceSBarry Smith .keywords: SNES, nonlinear, get, application, context 3749b94acceSBarry Smith 3759b94acceSBarry Smith .seealso: SNESSetApplicationContext() 3769b94acceSBarry Smith @*/ 3779b94acceSBarry Smith int SNESGetApplicationContext( SNES snes, void **usrP ) 3789b94acceSBarry Smith { 3793a40ed3dSBarry Smith PetscFunctionBegin; 38077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 3819b94acceSBarry Smith *usrP = snes->user; 3823a40ed3dSBarry Smith PetscFunctionReturn(0); 3839b94acceSBarry Smith } 38474679c65SBarry Smith 3855615d1e5SSatish Balay #undef __FUNC__ 386d4bb536fSBarry Smith #define __FUNC__ "SNESGetIterationNumber" 3879b94acceSBarry Smith /*@ 388c8228a4eSBarry Smith SNESGetIterationNumber - Gets the number of nonlinear iterations completed 389c8228a4eSBarry Smith at this time. 3909b94acceSBarry Smith 391c7afd0dbSLois Curfman McInnes Not Collective 392c7afd0dbSLois Curfman McInnes 3939b94acceSBarry Smith Input Parameter: 3949b94acceSBarry Smith . snes - SNES context 3959b94acceSBarry Smith 3969b94acceSBarry Smith Output Parameter: 3979b94acceSBarry Smith . iter - iteration number 3989b94acceSBarry Smith 399c8228a4eSBarry Smith Notes: 400c8228a4eSBarry Smith For example, during the computation of iteration 2 this would return 1. 401c8228a4eSBarry Smith 402c8228a4eSBarry Smith This is useful for using lagged Jacobians (where one does not recompute the 40308405cd6SLois Curfman McInnes Jacobian at each SNES iteration). For example, the code 40408405cd6SLois Curfman McInnes .vb 40508405cd6SLois Curfman McInnes ierr = SNESGetIterationNumber(snes,&it); 40608405cd6SLois Curfman McInnes if (!(it % 2)) { 40708405cd6SLois Curfman McInnes [compute Jacobian here] 40808405cd6SLois Curfman McInnes } 40908405cd6SLois Curfman McInnes .ve 410c8228a4eSBarry Smith can be used in your ComputeJacobian() function to cause the Jacobian to be 41108405cd6SLois Curfman McInnes recomputed every second SNES iteration. 412c8228a4eSBarry Smith 41336851e7fSLois Curfman McInnes Level: intermediate 41436851e7fSLois Curfman McInnes 4159b94acceSBarry Smith .keywords: SNES, nonlinear, get, iteration, number 4169b94acceSBarry Smith @*/ 4179b94acceSBarry Smith int SNESGetIterationNumber(SNES snes,int* iter) 4189b94acceSBarry Smith { 4193a40ed3dSBarry Smith PetscFunctionBegin; 42077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 42174679c65SBarry Smith PetscValidIntPointer(iter); 4229b94acceSBarry Smith *iter = snes->iter; 4233a40ed3dSBarry Smith PetscFunctionReturn(0); 4249b94acceSBarry Smith } 42574679c65SBarry Smith 4265615d1e5SSatish Balay #undef __FUNC__ 4275615d1e5SSatish Balay #define __FUNC__ "SNESGetFunctionNorm" 4289b94acceSBarry Smith /*@ 4299b94acceSBarry Smith SNESGetFunctionNorm - Gets the norm of the current function that was set 4309b94acceSBarry Smith with SNESSSetFunction(). 4319b94acceSBarry Smith 432c7afd0dbSLois Curfman McInnes Collective on SNES 433c7afd0dbSLois Curfman McInnes 4349b94acceSBarry Smith Input Parameter: 4359b94acceSBarry Smith . snes - SNES context 4369b94acceSBarry Smith 4379b94acceSBarry Smith Output Parameter: 4389b94acceSBarry Smith . fnorm - 2-norm of function 4399b94acceSBarry Smith 4409b94acceSBarry Smith Note: 4419b94acceSBarry Smith SNESGetFunctionNorm() is valid for SNES_NONLINEAR_EQUATIONS methods only. 442a86d99e1SLois Curfman McInnes A related routine for SNES_UNCONSTRAINED_MINIMIZATION methods is 443a86d99e1SLois Curfman McInnes SNESGetGradientNorm(). 4449b94acceSBarry Smith 44536851e7fSLois Curfman McInnes Level: intermediate 44636851e7fSLois Curfman McInnes 4479b94acceSBarry Smith .keywords: SNES, nonlinear, get, function, norm 448a86d99e1SLois Curfman McInnes 44908405cd6SLois Curfman McInnes .seealso: SNESGetFunction() 4509b94acceSBarry Smith @*/ 4519b94acceSBarry Smith int SNESGetFunctionNorm(SNES snes,Scalar *fnorm) 4529b94acceSBarry Smith { 4533a40ed3dSBarry Smith PetscFunctionBegin; 45477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 45574679c65SBarry Smith PetscValidScalarPointer(fnorm); 45674679c65SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 457d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_NONLINEAR_EQUATIONS only"); 45874679c65SBarry Smith } 4599b94acceSBarry Smith *fnorm = snes->norm; 4603a40ed3dSBarry Smith PetscFunctionReturn(0); 4619b94acceSBarry Smith } 46274679c65SBarry Smith 4635615d1e5SSatish Balay #undef __FUNC__ 4645615d1e5SSatish Balay #define __FUNC__ "SNESGetGradientNorm" 4659b94acceSBarry Smith /*@ 4669b94acceSBarry Smith SNESGetGradientNorm - Gets the norm of the current gradient that was set 4679b94acceSBarry Smith with SNESSSetGradient(). 4689b94acceSBarry Smith 469c7afd0dbSLois Curfman McInnes Collective on SNES 470c7afd0dbSLois Curfman McInnes 4719b94acceSBarry Smith Input Parameter: 4729b94acceSBarry Smith . snes - SNES context 4739b94acceSBarry Smith 4749b94acceSBarry Smith Output Parameter: 4759b94acceSBarry Smith . fnorm - 2-norm of gradient 4769b94acceSBarry Smith 4779b94acceSBarry Smith Note: 4789b94acceSBarry Smith SNESGetGradientNorm() is valid for SNES_UNCONSTRAINED_MINIMIZATION 479a86d99e1SLois Curfman McInnes methods only. A related routine for SNES_NONLINEAR_EQUATIONS methods 480a86d99e1SLois Curfman McInnes is SNESGetFunctionNorm(). 4819b94acceSBarry Smith 48236851e7fSLois Curfman McInnes Level: intermediate 48336851e7fSLois Curfman McInnes 4849b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient, norm 485a86d99e1SLois Curfman McInnes 486a86d99e1SLois Curfman McInnes .seelso: SNESSetGradient() 4879b94acceSBarry Smith @*/ 4889b94acceSBarry Smith int SNESGetGradientNorm(SNES snes,Scalar *gnorm) 4899b94acceSBarry Smith { 4903a40ed3dSBarry Smith PetscFunctionBegin; 49177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 49274679c65SBarry Smith PetscValidScalarPointer(gnorm); 49374679c65SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 494d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 49574679c65SBarry Smith } 4969b94acceSBarry Smith *gnorm = snes->norm; 4973a40ed3dSBarry Smith PetscFunctionReturn(0); 4989b94acceSBarry Smith } 49974679c65SBarry Smith 5005615d1e5SSatish Balay #undef __FUNC__ 501d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberUnsuccessfulSteps" 5029b94acceSBarry Smith /*@ 5039b94acceSBarry Smith SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 5049b94acceSBarry Smith attempted by the nonlinear solver. 5059b94acceSBarry Smith 506c7afd0dbSLois Curfman McInnes Not Collective 507c7afd0dbSLois Curfman McInnes 5089b94acceSBarry Smith Input Parameter: 5099b94acceSBarry Smith . snes - SNES context 5109b94acceSBarry Smith 5119b94acceSBarry Smith Output Parameter: 5129b94acceSBarry Smith . nfails - number of unsuccessful steps attempted 5139b94acceSBarry Smith 514c96a6f78SLois Curfman McInnes Notes: 515c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 516c96a6f78SLois Curfman McInnes 51736851e7fSLois Curfman McInnes Level: intermediate 51836851e7fSLois Curfman McInnes 5199b94acceSBarry Smith .keywords: SNES, nonlinear, get, number, unsuccessful, steps 5209b94acceSBarry Smith @*/ 5219b94acceSBarry Smith int SNESGetNumberUnsuccessfulSteps(SNES snes,int* nfails) 5229b94acceSBarry Smith { 5233a40ed3dSBarry Smith PetscFunctionBegin; 52477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 52574679c65SBarry Smith PetscValidIntPointer(nfails); 5269b94acceSBarry Smith *nfails = snes->nfailures; 5273a40ed3dSBarry Smith PetscFunctionReturn(0); 5289b94acceSBarry Smith } 529a847f771SSatish Balay 5305615d1e5SSatish Balay #undef __FUNC__ 531d4bb536fSBarry Smith #define __FUNC__ "SNESGetNumberLinearIterations" 532c96a6f78SLois Curfman McInnes /*@ 533c96a6f78SLois Curfman McInnes SNESGetNumberLinearIterations - Gets the total number of linear iterations 534c96a6f78SLois Curfman McInnes used by the nonlinear solver. 535c96a6f78SLois Curfman McInnes 536c7afd0dbSLois Curfman McInnes Not Collective 537c7afd0dbSLois Curfman McInnes 538c96a6f78SLois Curfman McInnes Input Parameter: 539c96a6f78SLois Curfman McInnes . snes - SNES context 540c96a6f78SLois Curfman McInnes 541c96a6f78SLois Curfman McInnes Output Parameter: 542c96a6f78SLois Curfman McInnes . lits - number of linear iterations 543c96a6f78SLois Curfman McInnes 544c96a6f78SLois Curfman McInnes Notes: 545c96a6f78SLois Curfman McInnes This counter is reset to zero for each successive call to SNESSolve(). 546c96a6f78SLois Curfman McInnes 54736851e7fSLois Curfman McInnes Level: intermediate 54836851e7fSLois Curfman McInnes 549c96a6f78SLois Curfman McInnes .keywords: SNES, nonlinear, get, number, linear, iterations 550c96a6f78SLois Curfman McInnes @*/ 55186bddb7dSBarry Smith int SNESGetNumberLinearIterations(SNES snes,int* lits) 552c96a6f78SLois Curfman McInnes { 5533a40ed3dSBarry Smith PetscFunctionBegin; 554c96a6f78SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 555c96a6f78SLois Curfman McInnes PetscValidIntPointer(lits); 556c96a6f78SLois Curfman McInnes *lits = snes->linear_its; 5573a40ed3dSBarry Smith PetscFunctionReturn(0); 558c96a6f78SLois Curfman McInnes } 559c96a6f78SLois Curfman McInnes 560c96a6f78SLois Curfman McInnes #undef __FUNC__ 561d4bb536fSBarry Smith #define __FUNC__ "SNESGetSLES" 5629b94acceSBarry Smith /*@C 5639b94acceSBarry Smith SNESGetSLES - Returns the SLES context for a SNES solver. 5649b94acceSBarry Smith 565c7afd0dbSLois Curfman McInnes Not Collective, but if SNES object is parallel, then SLES object is parallel 566c7afd0dbSLois Curfman McInnes 5679b94acceSBarry Smith Input Parameter: 5689b94acceSBarry Smith . snes - the SNES context 5699b94acceSBarry Smith 5709b94acceSBarry Smith Output Parameter: 5719b94acceSBarry Smith . sles - the SLES context 5729b94acceSBarry Smith 5739b94acceSBarry Smith Notes: 5749b94acceSBarry Smith The user can then directly manipulate the SLES context to set various 5759b94acceSBarry Smith options, etc. Likewise, the user can then extract and manipulate the 5769b94acceSBarry Smith KSP and PC contexts as well. 5779b94acceSBarry Smith 57836851e7fSLois Curfman McInnes Level: beginner 57936851e7fSLois Curfman McInnes 5809b94acceSBarry Smith .keywords: SNES, nonlinear, get, SLES, context 5819b94acceSBarry Smith 5829b94acceSBarry Smith .seealso: SLESGetPC(), SLESGetKSP() 5839b94acceSBarry Smith @*/ 5849b94acceSBarry Smith int SNESGetSLES(SNES snes,SLES *sles) 5859b94acceSBarry Smith { 5863a40ed3dSBarry Smith PetscFunctionBegin; 58777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 5889b94acceSBarry Smith *sles = snes->sles; 5893a40ed3dSBarry Smith PetscFunctionReturn(0); 5909b94acceSBarry Smith } 59182bf6240SBarry Smith 592e24b481bSBarry Smith #undef __FUNC__ 593e24b481bSBarry Smith #define __FUNC__ "SNESPublish_Petsc" 594*454a90a3SBarry Smith static int SNESPublish_Petsc(PetscObject obj) 595e24b481bSBarry Smith { 596aa482453SBarry Smith #if defined(PETSC_HAVE_AMS) 597*454a90a3SBarry Smith SNES v = (SNES) obj; 598e24b481bSBarry Smith int ierr; 59943d6d2cbSBarry Smith #endif 600e24b481bSBarry Smith 601e24b481bSBarry Smith PetscFunctionBegin; 602e24b481bSBarry Smith 60343d6d2cbSBarry Smith #if defined(PETSC_HAVE_AMS) 604e24b481bSBarry Smith /* if it is already published then return */ 605e24b481bSBarry Smith if (v->amem >=0 ) PetscFunctionReturn(0); 606e24b481bSBarry Smith 607*454a90a3SBarry Smith ierr = PetscObjectPublishBaseBegin(obj);CHKERRQ(ierr); 608e24b481bSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->iter,1,AMS_INT,AMS_READ, 609e24b481bSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 610e24b481bSBarry Smith ierr = AMS_Memory_add_field((AMS_Memory)v->amem,"Residual",&v->norm,1,AMS_DOUBLE,AMS_READ, 611e24b481bSBarry Smith AMS_COMMON,AMS_REDUCT_UNDEF);CHKERRQ(ierr); 612*454a90a3SBarry Smith ierr = PetscObjectPublishBaseEnd(obj);CHKERRQ(ierr); 613433994e6SBarry Smith #endif 614e24b481bSBarry Smith PetscFunctionReturn(0); 615e24b481bSBarry Smith } 616e24b481bSBarry Smith 6179b94acceSBarry Smith /* -----------------------------------------------------------*/ 6185615d1e5SSatish Balay #undef __FUNC__ 6195615d1e5SSatish Balay #define __FUNC__ "SNESCreate" 6209b94acceSBarry Smith /*@C 6219b94acceSBarry Smith SNESCreate - Creates a nonlinear solver context. 6229b94acceSBarry Smith 623c7afd0dbSLois Curfman McInnes Collective on MPI_Comm 624c7afd0dbSLois Curfman McInnes 625c7afd0dbSLois Curfman McInnes Input Parameters: 626c7afd0dbSLois Curfman McInnes + comm - MPI communicator 627c7afd0dbSLois Curfman McInnes - type - type of method, either 628c7afd0dbSLois Curfman McInnes SNES_NONLINEAR_EQUATIONS (for systems of nonlinear equations) 629c7afd0dbSLois Curfman McInnes or SNES_UNCONSTRAINED_MINIMIZATION (for unconstrained minimization) 6309b94acceSBarry Smith 6319b94acceSBarry Smith Output Parameter: 6329b94acceSBarry Smith . outsnes - the new SNES context 6339b94acceSBarry Smith 634c7afd0dbSLois Curfman McInnes Options Database Keys: 635c7afd0dbSLois Curfman McInnes + -snes_mf - Activates default matrix-free Jacobian-vector products, 636c7afd0dbSLois Curfman McInnes and no preconditioning matrix 637c7afd0dbSLois Curfman McInnes . -snes_mf_operator - Activates default matrix-free Jacobian-vector 638c7afd0dbSLois Curfman McInnes products, and a user-provided preconditioning matrix 639c7afd0dbSLois Curfman McInnes as set by SNESSetJacobian() 640c7afd0dbSLois Curfman McInnes - -snes_fd - Uses (slow!) finite differences to compute Jacobian 641c1f60f51SBarry Smith 64236851e7fSLois Curfman McInnes Level: beginner 64336851e7fSLois Curfman McInnes 6449b94acceSBarry Smith .keywords: SNES, nonlinear, create, context 6459b94acceSBarry Smith 64663a78c88SLois Curfman McInnes .seealso: SNESSolve(), SNESDestroy() 6479b94acceSBarry Smith @*/ 6484b0e389bSBarry Smith int SNESCreate(MPI_Comm comm,SNESProblemType type,SNES *outsnes) 6499b94acceSBarry Smith { 6509b94acceSBarry Smith int ierr; 6519b94acceSBarry Smith SNES snes; 6529b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 65337fcc0dbSBarry Smith 6543a40ed3dSBarry Smith PetscFunctionBegin; 6559b94acceSBarry Smith *outsnes = 0; 656d64ed03dSBarry Smith if (type != SNES_UNCONSTRAINED_MINIMIZATION && type != SNES_NONLINEAR_EQUATIONS){ 657d252947aSBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"incorrect method type"); 658d64ed03dSBarry Smith } 6593f1db9ecSBarry Smith PetscHeaderCreate(snes,_p_SNES,int,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView); 6609b94acceSBarry Smith PLogObjectCreate(snes); 661e24b481bSBarry Smith snes->bops->publish = SNESPublish_Petsc; 6629b94acceSBarry Smith snes->max_its = 50; 6639750a799SBarry Smith snes->max_funcs = 10000; 6649b94acceSBarry Smith snes->norm = 0.0; 6655a2d0531SBarry Smith if (type == SNES_UNCONSTRAINED_MINIMIZATION) { 666b18e04deSLois Curfman McInnes snes->rtol = 1.e-8; 667b18e04deSLois Curfman McInnes snes->ttol = 0.0; 6689b94acceSBarry Smith snes->atol = 1.e-10; 669d64ed03dSBarry Smith } else { 670b4874afaSBarry Smith snes->rtol = 1.e-8; 671b4874afaSBarry Smith snes->ttol = 0.0; 672b4874afaSBarry Smith snes->atol = 1.e-50; 673b4874afaSBarry Smith } 6749b94acceSBarry Smith snes->xtol = 1.e-8; 675b18e04deSLois Curfman McInnes snes->trunctol = 1.e-12; /* no longer used */ 6769b94acceSBarry Smith snes->nfuncs = 0; 6779b94acceSBarry Smith snes->nfailures = 0; 6787a00f4a9SLois Curfman McInnes snes->linear_its = 0; 679639f9d9dSBarry Smith snes->numbermonitors = 0; 6809b94acceSBarry Smith snes->data = 0; 6819b94acceSBarry Smith snes->view = 0; 6829b94acceSBarry Smith snes->computeumfunction = 0; 6839b94acceSBarry Smith snes->umfunP = 0; 6849b94acceSBarry Smith snes->fc = 0; 6859b94acceSBarry Smith snes->deltatol = 1.e-12; 6869b94acceSBarry Smith snes->fmin = -1.e30; 6879b94acceSBarry Smith snes->method_class = type; 6889b94acceSBarry Smith snes->set_method_called = 0; 68982bf6240SBarry Smith snes->setupcalled = 0; 6909b94acceSBarry Smith snes->ksp_ewconv = 0; 6916f24a144SLois Curfman McInnes snes->vwork = 0; 6926f24a144SLois Curfman McInnes snes->nwork = 0; 693758f92a0SBarry Smith snes->conv_hist_len = 0; 694758f92a0SBarry Smith snes->conv_hist_max = 0; 695758f92a0SBarry Smith snes->conv_hist = PETSC_NULL; 696758f92a0SBarry Smith snes->conv_hist_its = PETSC_NULL; 697758f92a0SBarry Smith snes->conv_hist_reset = PETSC_TRUE; 698184914b5SBarry Smith snes->reason = SNES_CONVERGED_ITERATING; 6999b94acceSBarry Smith 7009b94acceSBarry Smith /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 7010452661fSBarry Smith kctx = PetscNew(SNES_KSP_EW_ConvCtx);CHKPTRQ(kctx); 702eed86810SBarry Smith PLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 7039b94acceSBarry Smith snes->kspconvctx = (void*)kctx; 7049b94acceSBarry Smith kctx->version = 2; 7059b94acceSBarry Smith kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 7069b94acceSBarry Smith this was too large for some test cases */ 7079b94acceSBarry Smith kctx->rtol_last = 0; 7089b94acceSBarry Smith kctx->rtol_max = .9; 7099b94acceSBarry Smith kctx->gamma = 1.0; 7109b94acceSBarry Smith kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 7119b94acceSBarry Smith kctx->alpha = kctx->alpha2; 7129b94acceSBarry Smith kctx->threshold = .1; 7139b94acceSBarry Smith kctx->lresid_last = 0; 7149b94acceSBarry Smith kctx->norm_last = 0; 7159b94acceSBarry Smith 7169b94acceSBarry Smith ierr = SLESCreate(comm,&snes->sles);CHKERRQ(ierr); 7179b94acceSBarry Smith PLogObjectParent(snes,snes->sles) 7185334005bSBarry Smith 7199b94acceSBarry Smith *outsnes = snes; 720e24b481bSBarry Smith PetscPublishAll(snes); 7213a40ed3dSBarry Smith PetscFunctionReturn(0); 7229b94acceSBarry Smith } 7239b94acceSBarry Smith 7249b94acceSBarry Smith /* --------------------------------------------------------------- */ 7255615d1e5SSatish Balay #undef __FUNC__ 726d4bb536fSBarry Smith #define __FUNC__ "SNESSetFunction" 7279b94acceSBarry Smith /*@C 7289b94acceSBarry Smith SNESSetFunction - Sets the function evaluation routine and function 7299b94acceSBarry Smith vector for use by the SNES routines in solving systems of nonlinear 7309b94acceSBarry Smith equations. 7319b94acceSBarry Smith 732fee21e36SBarry Smith Collective on SNES 733fee21e36SBarry Smith 734c7afd0dbSLois Curfman McInnes Input Parameters: 735c7afd0dbSLois Curfman McInnes + snes - the SNES context 736c7afd0dbSLois Curfman McInnes . func - function evaluation routine 737c7afd0dbSLois Curfman McInnes . r - vector to store function value 738c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 739c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 7409b94acceSBarry Smith 741c7afd0dbSLois Curfman McInnes Calling sequence of func: 7428d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Vec f,void *ctx); 743c7afd0dbSLois Curfman McInnes 744313e4042SLois Curfman McInnes . f - function vector 745c7afd0dbSLois Curfman McInnes - ctx - optional user-defined function context 7469b94acceSBarry Smith 7479b94acceSBarry Smith Notes: 7489b94acceSBarry Smith The Newton-like methods typically solve linear systems of the form 7499b94acceSBarry Smith $ f'(x) x = -f(x), 750c7afd0dbSLois Curfman McInnes where f'(x) denotes the Jacobian matrix and f(x) is the function. 7519b94acceSBarry Smith 7529b94acceSBarry Smith SNESSetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 7539b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 7549b94acceSBarry Smith SNESSetMinimizationFunction() and SNESSetGradient(); 7559b94acceSBarry Smith 75636851e7fSLois Curfman McInnes Level: beginner 75736851e7fSLois Curfman McInnes 7589b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 7599b94acceSBarry Smith 760a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 7619b94acceSBarry Smith @*/ 7625334005bSBarry Smith int SNESSetFunction( SNES snes, Vec r, int (*func)(SNES,Vec,Vec,void*),void *ctx) 7639b94acceSBarry Smith { 7643a40ed3dSBarry Smith PetscFunctionBegin; 76577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 766184914b5SBarry Smith PetscValidHeaderSpecific(r,VEC_COOKIE); 767184914b5SBarry Smith PetscCheckSameComm(snes,r); 768a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 769a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 770a8c6a408SBarry Smith } 771184914b5SBarry Smith 7729b94acceSBarry Smith snes->computefunction = func; 7739b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 7749b94acceSBarry Smith snes->funP = ctx; 7753a40ed3dSBarry Smith PetscFunctionReturn(0); 7769b94acceSBarry Smith } 7779b94acceSBarry Smith 7785615d1e5SSatish Balay #undef __FUNC__ 7795615d1e5SSatish Balay #define __FUNC__ "SNESComputeFunction" 7809b94acceSBarry Smith /*@ 78136851e7fSLois Curfman McInnes SNESComputeFunction - Calls the function that has been set with 7829b94acceSBarry Smith SNESSetFunction(). 7839b94acceSBarry Smith 784c7afd0dbSLois Curfman McInnes Collective on SNES 785c7afd0dbSLois Curfman McInnes 7869b94acceSBarry Smith Input Parameters: 787c7afd0dbSLois Curfman McInnes + snes - the SNES context 788c7afd0dbSLois Curfman McInnes - x - input vector 7899b94acceSBarry Smith 7909b94acceSBarry Smith Output Parameter: 7913638b69dSLois Curfman McInnes . y - function vector, as set by SNESSetFunction() 7929b94acceSBarry Smith 7931bffabb2SLois Curfman McInnes Notes: 7949b94acceSBarry Smith SNESComputeFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only. 7959b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 7969b94acceSBarry Smith SNESComputeMinimizationFunction() and SNESComputeGradient(); 7979b94acceSBarry Smith 79836851e7fSLois Curfman McInnes SNESComputeFunction() is typically used within nonlinear solvers 79936851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 80036851e7fSLois Curfman McInnes themselves. 80136851e7fSLois Curfman McInnes 80236851e7fSLois Curfman McInnes Level: developer 80336851e7fSLois Curfman McInnes 8049b94acceSBarry Smith .keywords: SNES, nonlinear, compute, function 8059b94acceSBarry Smith 806a86d99e1SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetFunction() 8079b94acceSBarry Smith @*/ 8089b94acceSBarry Smith int SNESComputeFunction(SNES snes,Vec x, Vec y) 8099b94acceSBarry Smith { 8109b94acceSBarry Smith int ierr; 8119b94acceSBarry Smith 8123a40ed3dSBarry Smith PetscFunctionBegin; 813184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 814184914b5SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 815184914b5SBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE); 816184914b5SBarry Smith PetscCheckSameComm(snes,x); 817184914b5SBarry Smith PetscCheckSameComm(snes,y); 818d4bb536fSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 819a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 820d4bb536fSBarry Smith } 821184914b5SBarry Smith 8229b94acceSBarry Smith PLogEventBegin(SNES_FunctionEval,snes,x,y,0); 823d64ed03dSBarry Smith PetscStackPush("SNES user function"); 8249b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 825d64ed03dSBarry Smith PetscStackPop; 826ae3c334cSLois Curfman McInnes snes->nfuncs++; 8279b94acceSBarry Smith PLogEventEnd(SNES_FunctionEval,snes,x,y,0); 8283a40ed3dSBarry Smith PetscFunctionReturn(0); 8299b94acceSBarry Smith } 8309b94acceSBarry Smith 8315615d1e5SSatish Balay #undef __FUNC__ 832d4bb536fSBarry Smith #define __FUNC__ "SNESSetMinimizationFunction" 8339b94acceSBarry Smith /*@C 8349b94acceSBarry Smith SNESSetMinimizationFunction - Sets the function evaluation routine for 8359b94acceSBarry Smith unconstrained minimization. 8369b94acceSBarry Smith 837fee21e36SBarry Smith Collective on SNES 838fee21e36SBarry Smith 839c7afd0dbSLois Curfman McInnes Input Parameters: 840c7afd0dbSLois Curfman McInnes + snes - the SNES context 841c7afd0dbSLois Curfman McInnes . func - function evaluation routine 842c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 843c7afd0dbSLois Curfman McInnes function evaluation routine (may be PETSC_NULL) 8449b94acceSBarry Smith 845c7afd0dbSLois Curfman McInnes Calling sequence of func: 8468d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,double *f,void *ctx); 847c7afd0dbSLois Curfman McInnes 848c7afd0dbSLois Curfman McInnes + x - input vector 8499b94acceSBarry Smith . f - function 850c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined function context 8519b94acceSBarry Smith 85236851e7fSLois Curfman McInnes Level: beginner 85336851e7fSLois Curfman McInnes 8549b94acceSBarry Smith Notes: 8559b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 8569b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 8579b94acceSBarry Smith SNESSetFunction(). 8589b94acceSBarry Smith 8599b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimization, function 8609b94acceSBarry Smith 861a86d99e1SLois Curfman McInnes .seealso: SNESGetMinimizationFunction(), SNESComputeMinimizationFunction(), 862a86d99e1SLois Curfman McInnes SNESSetHessian(), SNESSetGradient() 8639b94acceSBarry Smith @*/ 864*454a90a3SBarry Smith int SNESSetMinimizationFunction(SNES snes,int (*func)(SNES,Vec,double*,void*),void *ctx) 8659b94acceSBarry Smith { 8663a40ed3dSBarry Smith PetscFunctionBegin; 86777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 868a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 869a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 870a8c6a408SBarry Smith } 8719b94acceSBarry Smith snes->computeumfunction = func; 8729b94acceSBarry Smith snes->umfunP = ctx; 8733a40ed3dSBarry Smith PetscFunctionReturn(0); 8749b94acceSBarry Smith } 8759b94acceSBarry Smith 8765615d1e5SSatish Balay #undef __FUNC__ 8775615d1e5SSatish Balay #define __FUNC__ "SNESComputeMinimizationFunction" 8789b94acceSBarry Smith /*@ 8799b94acceSBarry Smith SNESComputeMinimizationFunction - Computes the function that has been 8809b94acceSBarry Smith set with SNESSetMinimizationFunction(). 8819b94acceSBarry Smith 882c7afd0dbSLois Curfman McInnes Collective on SNES 883c7afd0dbSLois Curfman McInnes 8849b94acceSBarry Smith Input Parameters: 885c7afd0dbSLois Curfman McInnes + snes - the SNES context 886c7afd0dbSLois Curfman McInnes - x - input vector 8879b94acceSBarry Smith 8889b94acceSBarry Smith Output Parameter: 8899b94acceSBarry Smith . y - function value 8909b94acceSBarry Smith 8919b94acceSBarry Smith Notes: 8929b94acceSBarry Smith SNESComputeMinimizationFunction() is valid only for 8939b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 8949b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 895a86d99e1SLois Curfman McInnes 89636851e7fSLois Curfman McInnes SNESComputeMinimizationFunction() is typically used within minimization 89736851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 89836851e7fSLois Curfman McInnes themselves. 89936851e7fSLois Curfman McInnes 90036851e7fSLois Curfman McInnes Level: developer 90136851e7fSLois Curfman McInnes 902a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, minimization, function 903a86d99e1SLois Curfman McInnes 904a86d99e1SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESGetMinimizationFunction(), 905a86d99e1SLois Curfman McInnes SNESComputeGradient(), SNESComputeHessian() 9069b94acceSBarry Smith @*/ 9079b94acceSBarry Smith int SNESComputeMinimizationFunction(SNES snes,Vec x,double *y) 9089b94acceSBarry Smith { 9099b94acceSBarry Smith int ierr; 9103a40ed3dSBarry Smith 9113a40ed3dSBarry Smith PetscFunctionBegin; 912184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 913184914b5SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 914184914b5SBarry Smith PetscCheckSameComm(snes,x); 915a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 916a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"Only for SNES_UNCONSTRAINED_MINIMIZATION"); 917a8c6a408SBarry Smith } 918184914b5SBarry Smith 9199b94acceSBarry Smith PLogEventBegin(SNES_MinimizationFunctionEval,snes,x,y,0); 920d64ed03dSBarry Smith PetscStackPush("SNES user minimzation function"); 9219b94acceSBarry Smith ierr = (*snes->computeumfunction)(snes,x,y,snes->umfunP);CHKERRQ(ierr); 922d64ed03dSBarry Smith PetscStackPop; 923ae3c334cSLois Curfman McInnes snes->nfuncs++; 9249b94acceSBarry Smith PLogEventEnd(SNES_MinimizationFunctionEval,snes,x,y,0); 9253a40ed3dSBarry Smith PetscFunctionReturn(0); 9269b94acceSBarry Smith } 9279b94acceSBarry Smith 9285615d1e5SSatish Balay #undef __FUNC__ 929d4bb536fSBarry Smith #define __FUNC__ "SNESSetGradient" 9309b94acceSBarry Smith /*@C 9319b94acceSBarry Smith SNESSetGradient - Sets the gradient evaluation routine and gradient 9329b94acceSBarry Smith vector for use by the SNES routines. 9339b94acceSBarry Smith 934c7afd0dbSLois Curfman McInnes Collective on SNES 935c7afd0dbSLois Curfman McInnes 9369b94acceSBarry Smith Input Parameters: 937c7afd0dbSLois Curfman McInnes + snes - the SNES context 9389b94acceSBarry Smith . func - function evaluation routine 939044dda88SLois Curfman McInnes . ctx - optional user-defined context for private data for the 940044dda88SLois Curfman McInnes gradient evaluation routine (may be PETSC_NULL) 941c7afd0dbSLois Curfman McInnes - r - vector to store gradient value 942fee21e36SBarry Smith 9439b94acceSBarry Smith Calling sequence of func: 9448d76a1e5SLois Curfman McInnes $ func (SNES, Vec x, Vec g, void *ctx); 9459b94acceSBarry Smith 946c7afd0dbSLois Curfman McInnes + x - input vector 9479b94acceSBarry Smith . g - gradient vector 948c7afd0dbSLois Curfman McInnes - ctx - optional user-defined gradient context 9499b94acceSBarry Smith 9509b94acceSBarry Smith Notes: 9519b94acceSBarry Smith SNESSetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 9529b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 9539b94acceSBarry Smith SNESSetFunction(). 9549b94acceSBarry Smith 95536851e7fSLois Curfman McInnes Level: beginner 95636851e7fSLois Curfman McInnes 9579b94acceSBarry Smith .keywords: SNES, nonlinear, set, function 9589b94acceSBarry Smith 959a86d99e1SLois Curfman McInnes .seealso: SNESGetGradient(), SNESComputeGradient(), SNESSetHessian(), 960a86d99e1SLois Curfman McInnes SNESSetMinimizationFunction(), 9619b94acceSBarry Smith @*/ 96274679c65SBarry Smith int SNESSetGradient(SNES snes,Vec r,int (*func)(SNES,Vec,Vec,void*),void *ctx) 9639b94acceSBarry Smith { 9643a40ed3dSBarry Smith PetscFunctionBegin; 96577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 966184914b5SBarry Smith PetscValidHeaderSpecific(r,VEC_COOKIE); 967184914b5SBarry Smith PetscCheckSameComm(snes,r); 968a8c6a408SBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 969a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 970a8c6a408SBarry Smith } 9719b94acceSBarry Smith snes->computefunction = func; 9729b94acceSBarry Smith snes->vec_func = snes->vec_func_always = r; 9739b94acceSBarry Smith snes->funP = ctx; 9743a40ed3dSBarry Smith PetscFunctionReturn(0); 9759b94acceSBarry Smith } 9769b94acceSBarry Smith 9775615d1e5SSatish Balay #undef __FUNC__ 9785615d1e5SSatish Balay #define __FUNC__ "SNESComputeGradient" 9799b94acceSBarry Smith /*@ 980a86d99e1SLois Curfman McInnes SNESComputeGradient - Computes the gradient that has been set with 981a86d99e1SLois Curfman McInnes SNESSetGradient(). 9829b94acceSBarry Smith 983c7afd0dbSLois Curfman McInnes Collective on SNES 984c7afd0dbSLois Curfman McInnes 9859b94acceSBarry Smith Input Parameters: 986c7afd0dbSLois Curfman McInnes + snes - the SNES context 987c7afd0dbSLois Curfman McInnes - x - input vector 9889b94acceSBarry Smith 9899b94acceSBarry Smith Output Parameter: 9909b94acceSBarry Smith . y - gradient vector 9919b94acceSBarry Smith 9929b94acceSBarry Smith Notes: 9939b94acceSBarry Smith SNESComputeGradient() is valid only for 9949b94acceSBarry Smith SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 9959b94acceSBarry Smith SNES_NONLINEAR_EQUATIONS methods is SNESComputeFunction(). 996a86d99e1SLois Curfman McInnes 99736851e7fSLois Curfman McInnes SNESComputeGradient() is typically used within minimization 99836851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 99936851e7fSLois Curfman McInnes themselves. 100036851e7fSLois Curfman McInnes 100136851e7fSLois Curfman McInnes Level: developer 100236851e7fSLois Curfman McInnes 1003a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, compute, gradient 1004a86d99e1SLois Curfman McInnes 1005a86d99e1SLois Curfman McInnes .seealso: SNESSetGradient(), SNESGetGradient(), 1006a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction(), SNESComputeHessian() 10079b94acceSBarry Smith @*/ 10089b94acceSBarry Smith int SNESComputeGradient(SNES snes,Vec x, Vec y) 10099b94acceSBarry Smith { 10109b94acceSBarry Smith int ierr; 10113a40ed3dSBarry Smith 10123a40ed3dSBarry Smith PetscFunctionBegin; 1013184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1014184914b5SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 1015184914b5SBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE); 1016184914b5SBarry Smith PetscCheckSameComm(snes,x); 1017184914b5SBarry Smith PetscCheckSameComm(snes,y); 10183a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1019a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 10203a40ed3dSBarry Smith } 10219b94acceSBarry Smith PLogEventBegin(SNES_GradientEval,snes,x,y,0); 1022d64ed03dSBarry Smith PetscStackPush("SNES user gradient function"); 10239b94acceSBarry Smith ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 1024d64ed03dSBarry Smith PetscStackPop; 10259b94acceSBarry Smith PLogEventEnd(SNES_GradientEval,snes,x,y,0); 10263a40ed3dSBarry Smith PetscFunctionReturn(0); 10279b94acceSBarry Smith } 10289b94acceSBarry Smith 10295615d1e5SSatish Balay #undef __FUNC__ 10305615d1e5SSatish Balay #define __FUNC__ "SNESComputeJacobian" 103162fef451SLois Curfman McInnes /*@ 103262fef451SLois Curfman McInnes SNESComputeJacobian - Computes the Jacobian matrix that has been 103362fef451SLois Curfman McInnes set with SNESSetJacobian(). 103462fef451SLois Curfman McInnes 1035c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1036c7afd0dbSLois Curfman McInnes 103762fef451SLois Curfman McInnes Input Parameters: 1038c7afd0dbSLois Curfman McInnes + snes - the SNES context 1039c7afd0dbSLois Curfman McInnes - x - input vector 104062fef451SLois Curfman McInnes 104162fef451SLois Curfman McInnes Output Parameters: 1042c7afd0dbSLois Curfman McInnes + A - Jacobian matrix 104362fef451SLois Curfman McInnes . B - optional preconditioning matrix 1044c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 1045fee21e36SBarry Smith 104662fef451SLois Curfman McInnes Notes: 104762fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 104862fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 104962fef451SLois Curfman McInnes 1050dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 1051dc5a77f8SLois Curfman McInnes flag parameter. 105262fef451SLois Curfman McInnes 105362fef451SLois Curfman McInnes SNESComputeJacobian() is valid only for SNES_NONLINEAR_EQUATIONS 105462fef451SLois Curfman McInnes methods. An analogous routine for SNES_UNCONSTRAINED_MINIMIZATION 1055005c665bSBarry Smith methods is SNESComputeHessian(). 105662fef451SLois Curfman McInnes 105736851e7fSLois Curfman McInnes SNESComputeJacobian() is typically used within nonlinear solver 105836851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 105936851e7fSLois Curfman McInnes themselves. 106036851e7fSLois Curfman McInnes 106136851e7fSLois Curfman McInnes Level: developer 106236851e7fSLois Curfman McInnes 106362fef451SLois Curfman McInnes .keywords: SNES, compute, Jacobian, matrix 106462fef451SLois Curfman McInnes 106562fef451SLois Curfman McInnes .seealso: SNESSetJacobian(), SLESSetOperators() 106662fef451SLois Curfman McInnes @*/ 10679b94acceSBarry Smith int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 10689b94acceSBarry Smith { 10699b94acceSBarry Smith int ierr; 10703a40ed3dSBarry Smith 10713a40ed3dSBarry Smith PetscFunctionBegin; 1072184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1073184914b5SBarry Smith PetscValidHeaderSpecific(X,VEC_COOKIE); 1074184914b5SBarry Smith PetscCheckSameComm(snes,X); 10753a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1076a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 10773a40ed3dSBarry Smith } 10783a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 10799b94acceSBarry Smith PLogEventBegin(SNES_JacobianEval,snes,X,*A,*B); 1080c4fc05e7SBarry Smith *flg = DIFFERENT_NONZERO_PATTERN; 1081d64ed03dSBarry Smith PetscStackPush("SNES user Jacobian function"); 10829b94acceSBarry Smith ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 1083d64ed03dSBarry Smith PetscStackPop; 10849b94acceSBarry Smith PLogEventEnd(SNES_JacobianEval,snes,X,*A,*B); 10856d84be18SBarry Smith /* make sure user returned a correct Jacobian and preconditioner */ 108677c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 108777c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 10883a40ed3dSBarry Smith PetscFunctionReturn(0); 10899b94acceSBarry Smith } 10909b94acceSBarry Smith 10915615d1e5SSatish Balay #undef __FUNC__ 10925615d1e5SSatish Balay #define __FUNC__ "SNESComputeHessian" 109362fef451SLois Curfman McInnes /*@ 109462fef451SLois Curfman McInnes SNESComputeHessian - Computes the Hessian matrix that has been 109562fef451SLois Curfman McInnes set with SNESSetHessian(). 109662fef451SLois Curfman McInnes 1097c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1098c7afd0dbSLois Curfman McInnes 109962fef451SLois Curfman McInnes Input Parameters: 1100c7afd0dbSLois Curfman McInnes + snes - the SNES context 1101c7afd0dbSLois Curfman McInnes - x - input vector 110262fef451SLois Curfman McInnes 110362fef451SLois Curfman McInnes Output Parameters: 1104c7afd0dbSLois Curfman McInnes + A - Hessian matrix 110562fef451SLois Curfman McInnes . B - optional preconditioning matrix 1106c7afd0dbSLois Curfman McInnes - flag - flag indicating matrix structure 1107fee21e36SBarry Smith 110862fef451SLois Curfman McInnes Notes: 110962fef451SLois Curfman McInnes Most users should not need to explicitly call this routine, as it 111062fef451SLois Curfman McInnes is used internally within the nonlinear solvers. 111162fef451SLois Curfman McInnes 1112dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the 1113dc5a77f8SLois Curfman McInnes flag parameter. 111462fef451SLois Curfman McInnes 111562fef451SLois Curfman McInnes SNESComputeHessian() is valid only for 111662fef451SLois Curfman McInnes SNES_UNCONSTRAINED_MINIMIZATION methods. An analogous routine for 111762fef451SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods is SNESComputeJacobian(). 111862fef451SLois Curfman McInnes 111936851e7fSLois Curfman McInnes SNESComputeHessian() is typically used within minimization 112036851e7fSLois Curfman McInnes implementations, so most users would not generally call this routine 112136851e7fSLois Curfman McInnes themselves. 112236851e7fSLois Curfman McInnes 112336851e7fSLois Curfman McInnes Level: developer 112436851e7fSLois Curfman McInnes 112562fef451SLois Curfman McInnes .keywords: SNES, compute, Hessian, matrix 112662fef451SLois Curfman McInnes 1127a86d99e1SLois Curfman McInnes .seealso: SNESSetHessian(), SLESSetOperators(), SNESComputeGradient(), 1128a86d99e1SLois Curfman McInnes SNESComputeMinimizationFunction() 112962fef451SLois Curfman McInnes @*/ 113062fef451SLois Curfman McInnes int SNESComputeHessian(SNES snes,Vec x,Mat *A,Mat *B,MatStructure *flag) 11319b94acceSBarry Smith { 11329b94acceSBarry Smith int ierr; 11333a40ed3dSBarry Smith 11343a40ed3dSBarry Smith PetscFunctionBegin; 1135184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1136184914b5SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 1137184914b5SBarry Smith PetscCheckSameComm(snes,x); 11383a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1139a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 11403a40ed3dSBarry Smith } 11413a40ed3dSBarry Smith if (!snes->computejacobian) PetscFunctionReturn(0); 114262fef451SLois Curfman McInnes PLogEventBegin(SNES_HessianEval,snes,x,*A,*B); 11430f4a323eSLois Curfman McInnes *flag = DIFFERENT_NONZERO_PATTERN; 1144d64ed03dSBarry Smith PetscStackPush("SNES user Hessian function"); 114562fef451SLois Curfman McInnes ierr = (*snes->computejacobian)(snes,x,A,B,flag,snes->jacP);CHKERRQ(ierr); 1146d64ed03dSBarry Smith PetscStackPop; 114762fef451SLois Curfman McInnes PLogEventEnd(SNES_HessianEval,snes,x,*A,*B); 11480f4a323eSLois Curfman McInnes /* make sure user returned a correct Jacobian and preconditioner */ 114977c4ece6SBarry Smith PetscValidHeaderSpecific(*A,MAT_COOKIE); 115077c4ece6SBarry Smith PetscValidHeaderSpecific(*B,MAT_COOKIE); 11513a40ed3dSBarry Smith PetscFunctionReturn(0); 11529b94acceSBarry Smith } 11539b94acceSBarry Smith 11545615d1e5SSatish Balay #undef __FUNC__ 1155d4bb536fSBarry Smith #define __FUNC__ "SNESSetJacobian" 11569b94acceSBarry Smith /*@C 11579b94acceSBarry Smith SNESSetJacobian - Sets the function to compute Jacobian as well as the 1158044dda88SLois Curfman McInnes location to store the matrix. 11599b94acceSBarry Smith 1160c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1161c7afd0dbSLois Curfman McInnes 11629b94acceSBarry Smith Input Parameters: 1163c7afd0dbSLois Curfman McInnes + snes - the SNES context 11649b94acceSBarry Smith . A - Jacobian matrix 11659b94acceSBarry Smith . B - preconditioner matrix (usually same as the Jacobian) 11669b94acceSBarry Smith . func - Jacobian evaluation routine 1167c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 11682cd2dfdcSLois Curfman McInnes Jacobian evaluation routine (may be PETSC_NULL) 11699b94acceSBarry Smith 11709b94acceSBarry Smith Calling sequence of func: 11718d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 11729b94acceSBarry Smith 1173c7afd0dbSLois Curfman McInnes + x - input vector 11749b94acceSBarry Smith . A - Jacobian matrix 11759b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1176ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1177ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 1178c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Jacobian context 11799b94acceSBarry Smith 11809b94acceSBarry Smith Notes: 1181dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 11822cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1183ac21db08SLois Curfman McInnes 1184ac21db08SLois Curfman McInnes The routine func() takes Mat * as the matrix arguments rather than Mat. 11859b94acceSBarry Smith This allows the Jacobian evaluation routine to replace A and/or B with a 11869b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 11879b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 11889b94acceSBarry Smith throughout the global iterations. 11899b94acceSBarry Smith 119036851e7fSLois Curfman McInnes Level: beginner 119136851e7fSLois Curfman McInnes 11929b94acceSBarry Smith .keywords: SNES, nonlinear, set, Jacobian, matrix 11939b94acceSBarry Smith 1194ac21db08SLois Curfman McInnes .seealso: SLESSetOperators(), SNESSetFunction() 11959b94acceSBarry Smith @*/ 1196*454a90a3SBarry Smith int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 11979b94acceSBarry Smith { 11983a40ed3dSBarry Smith PetscFunctionBegin; 119977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1200184914b5SBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE); 1201184914b5SBarry Smith PetscValidHeaderSpecific(B,MAT_COOKIE); 1202184914b5SBarry Smith PetscCheckSameComm(snes,A); 1203184914b5SBarry Smith PetscCheckSameComm(snes,B); 1204a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1205a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 1206a8c6a408SBarry Smith } 1207184914b5SBarry Smith 12089b94acceSBarry Smith snes->computejacobian = func; 12099b94acceSBarry Smith snes->jacP = ctx; 12109b94acceSBarry Smith snes->jacobian = A; 12119b94acceSBarry Smith snes->jacobian_pre = B; 12123a40ed3dSBarry Smith PetscFunctionReturn(0); 12139b94acceSBarry Smith } 121462fef451SLois Curfman McInnes 12155615d1e5SSatish Balay #undef __FUNC__ 1216d4bb536fSBarry Smith #define __FUNC__ "SNESGetJacobian" 1217c2aafc4cSSatish Balay /*@C 1218b4fd4287SBarry Smith SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1219b4fd4287SBarry Smith provided context for evaluating the Jacobian. 1220b4fd4287SBarry Smith 1221c7afd0dbSLois Curfman McInnes Not Collective, but Mat object will be parallel if SNES object is 1222c7afd0dbSLois Curfman McInnes 1223b4fd4287SBarry Smith Input Parameter: 1224b4fd4287SBarry Smith . snes - the nonlinear solver context 1225b4fd4287SBarry Smith 1226b4fd4287SBarry Smith Output Parameters: 1227c7afd0dbSLois Curfman McInnes + A - location to stash Jacobian matrix (or PETSC_NULL) 1228b4fd4287SBarry Smith . B - location to stash preconditioner matrix (or PETSC_NULL) 1229c7afd0dbSLois Curfman McInnes - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1230fee21e36SBarry Smith 123136851e7fSLois Curfman McInnes Level: advanced 123236851e7fSLois Curfman McInnes 1233b4fd4287SBarry Smith .seealso: SNESSetJacobian(), SNESComputeJacobian() 1234b4fd4287SBarry Smith @*/ 1235b4fd4287SBarry Smith int SNESGetJacobian(SNES snes,Mat *A,Mat *B, void **ctx) 1236b4fd4287SBarry Smith { 12373a40ed3dSBarry Smith PetscFunctionBegin; 1238184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 12393a40ed3dSBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 1240a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 12413a40ed3dSBarry Smith } 1242b4fd4287SBarry Smith if (A) *A = snes->jacobian; 1243b4fd4287SBarry Smith if (B) *B = snes->jacobian_pre; 1244b4fd4287SBarry Smith if (ctx) *ctx = snes->jacP; 12453a40ed3dSBarry Smith PetscFunctionReturn(0); 1246b4fd4287SBarry Smith } 1247b4fd4287SBarry Smith 12485615d1e5SSatish Balay #undef __FUNC__ 1249d4bb536fSBarry Smith #define __FUNC__ "SNESSetHessian" 12509b94acceSBarry Smith /*@C 12519b94acceSBarry Smith SNESSetHessian - Sets the function to compute Hessian as well as the 1252044dda88SLois Curfman McInnes location to store the matrix. 12539b94acceSBarry Smith 1254c7afd0dbSLois Curfman McInnes Collective on SNES and Mat 1255c7afd0dbSLois Curfman McInnes 12569b94acceSBarry Smith Input Parameters: 1257c7afd0dbSLois Curfman McInnes + snes - the SNES context 12589b94acceSBarry Smith . A - Hessian matrix 12599b94acceSBarry Smith . B - preconditioner matrix (usually same as the Hessian) 12609b94acceSBarry Smith . func - Jacobian evaluation routine 1261c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined context for private data for the 1262313e4042SLois Curfman McInnes Hessian evaluation routine (may be PETSC_NULL) 12639b94acceSBarry Smith 12649b94acceSBarry Smith Calling sequence of func: 12658d76a1e5SLois Curfman McInnes $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 12669b94acceSBarry Smith 1267c7afd0dbSLois Curfman McInnes + x - input vector 12689b94acceSBarry Smith . A - Hessian matrix 12699b94acceSBarry Smith . B - preconditioner matrix, usually the same as A 1270ac21db08SLois Curfman McInnes . flag - flag indicating information about the preconditioner matrix 1271ac21db08SLois Curfman McInnes structure (same as flag in SLESSetOperators()) 1272c7afd0dbSLois Curfman McInnes - ctx - [optional] user-defined Hessian context 12739b94acceSBarry Smith 12749b94acceSBarry Smith Notes: 1275dc5a77f8SLois Curfman McInnes See SLESSetOperators() for important information about setting the flag 12762cd2dfdcSLois Curfman McInnes output parameter in the routine func(). Be sure to read this information! 1277ac21db08SLois Curfman McInnes 12789b94acceSBarry Smith The function func() takes Mat * as the matrix arguments rather than Mat. 12799b94acceSBarry Smith This allows the Hessian evaluation routine to replace A and/or B with a 12809b94acceSBarry Smith completely new new matrix structure (not just different matrix elements) 12819b94acceSBarry Smith when appropriate, for instance, if the nonzero structure is changing 12829b94acceSBarry Smith throughout the global iterations. 12839b94acceSBarry Smith 128436851e7fSLois Curfman McInnes Level: beginner 128536851e7fSLois Curfman McInnes 12869b94acceSBarry Smith .keywords: SNES, nonlinear, set, Hessian, matrix 12879b94acceSBarry Smith 1288ac21db08SLois Curfman McInnes .seealso: SNESSetMinimizationFunction(), SNESSetGradient(), SLESSetOperators() 12899b94acceSBarry Smith @*/ 1290*454a90a3SBarry Smith int SNESSetHessian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 12919b94acceSBarry Smith { 12923a40ed3dSBarry Smith PetscFunctionBegin; 129377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1294184914b5SBarry Smith PetscValidHeaderSpecific(A,MAT_COOKIE); 1295184914b5SBarry Smith PetscValidHeaderSpecific(B,MAT_COOKIE); 1296184914b5SBarry Smith PetscCheckSameComm(snes,A); 1297184914b5SBarry Smith PetscCheckSameComm(snes,B); 1298d4bb536fSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 1299a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 1300d4bb536fSBarry Smith } 13019b94acceSBarry Smith snes->computejacobian = func; 13029b94acceSBarry Smith snes->jacP = ctx; 13039b94acceSBarry Smith snes->jacobian = A; 13049b94acceSBarry Smith snes->jacobian_pre = B; 13053a40ed3dSBarry Smith PetscFunctionReturn(0); 13069b94acceSBarry Smith } 13079b94acceSBarry Smith 13085615d1e5SSatish Balay #undef __FUNC__ 1309d4bb536fSBarry Smith #define __FUNC__ "SNESGetHessian" 131062fef451SLois Curfman McInnes /*@ 131162fef451SLois Curfman McInnes SNESGetHessian - Returns the Hessian matrix and optionally the user 131262fef451SLois Curfman McInnes provided context for evaluating the Hessian. 131362fef451SLois Curfman McInnes 1314c7afd0dbSLois Curfman McInnes Not Collective, but Mat object is parallel if SNES object is parallel 1315c7afd0dbSLois Curfman McInnes 131662fef451SLois Curfman McInnes Input Parameter: 131762fef451SLois Curfman McInnes . snes - the nonlinear solver context 131862fef451SLois Curfman McInnes 131962fef451SLois Curfman McInnes Output Parameters: 1320c7afd0dbSLois Curfman McInnes + A - location to stash Hessian matrix (or PETSC_NULL) 132162fef451SLois Curfman McInnes . B - location to stash preconditioner matrix (or PETSC_NULL) 1322c7afd0dbSLois Curfman McInnes - ctx - location to stash Hessian ctx (or PETSC_NULL) 1323fee21e36SBarry Smith 132436851e7fSLois Curfman McInnes Level: advanced 132536851e7fSLois Curfman McInnes 132662fef451SLois Curfman McInnes .seealso: SNESSetHessian(), SNESComputeHessian() 1327c7afd0dbSLois Curfman McInnes 1328c7afd0dbSLois Curfman McInnes .keywords: SNES, get, Hessian 132962fef451SLois Curfman McInnes @*/ 133062fef451SLois Curfman McInnes int SNESGetHessian(SNES snes,Mat *A,Mat *B, void **ctx) 133162fef451SLois Curfman McInnes { 13323a40ed3dSBarry Smith PetscFunctionBegin; 1333184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 13343a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION){ 1335a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 13363a40ed3dSBarry Smith } 133762fef451SLois Curfman McInnes if (A) *A = snes->jacobian; 133862fef451SLois Curfman McInnes if (B) *B = snes->jacobian_pre; 133962fef451SLois Curfman McInnes if (ctx) *ctx = snes->jacP; 13403a40ed3dSBarry Smith PetscFunctionReturn(0); 134162fef451SLois Curfman McInnes } 134262fef451SLois Curfman McInnes 13439b94acceSBarry Smith /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 13449b94acceSBarry Smith 13455615d1e5SSatish Balay #undef __FUNC__ 13465615d1e5SSatish Balay #define __FUNC__ "SNESSetUp" 13479b94acceSBarry Smith /*@ 13489b94acceSBarry Smith SNESSetUp - Sets up the internal data structures for the later use 1349272ac6f2SLois Curfman McInnes of a nonlinear solver. 13509b94acceSBarry Smith 1351fee21e36SBarry Smith Collective on SNES 1352fee21e36SBarry Smith 1353c7afd0dbSLois Curfman McInnes Input Parameters: 1354c7afd0dbSLois Curfman McInnes + snes - the SNES context 1355c7afd0dbSLois Curfman McInnes - x - the solution vector 1356c7afd0dbSLois Curfman McInnes 1357272ac6f2SLois Curfman McInnes Notes: 1358272ac6f2SLois Curfman McInnes For basic use of the SNES solvers the user need not explicitly call 1359272ac6f2SLois Curfman McInnes SNESSetUp(), since these actions will automatically occur during 1360272ac6f2SLois Curfman McInnes the call to SNESSolve(). However, if one wishes to control this 1361272ac6f2SLois Curfman McInnes phase separately, SNESSetUp() should be called after SNESCreate() 1362272ac6f2SLois Curfman McInnes and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1363272ac6f2SLois Curfman McInnes 136436851e7fSLois Curfman McInnes Level: advanced 136536851e7fSLois Curfman McInnes 13669b94acceSBarry Smith .keywords: SNES, nonlinear, setup 13679b94acceSBarry Smith 13689b94acceSBarry Smith .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 13699b94acceSBarry Smith @*/ 13708ddd3da0SLois Curfman McInnes int SNESSetUp(SNES snes,Vec x) 13719b94acceSBarry Smith { 1372272ac6f2SLois Curfman McInnes int ierr, flg; 13733a40ed3dSBarry Smith 13743a40ed3dSBarry Smith PetscFunctionBegin; 137577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 137677c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 1377184914b5SBarry Smith PetscCheckSameComm(snes,x); 13788ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 13799b94acceSBarry Smith 1380c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);CHKERRQ(ierr); 1381c1f60f51SBarry Smith /* 1382c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1383dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1384c1f60f51SBarry Smith */ 1385c1f60f51SBarry Smith if (flg) { 1386c1f60f51SBarry Smith Mat J; 13875a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1388c1f60f51SBarry Smith PLogObjectParent(snes,J); 1389c1f60f51SBarry Smith snes->mfshell = J; 1390c1f60f51SBarry Smith snes->jacobian = J; 1391c2aafc4cSSatish Balay if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1392c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1393d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1394c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1395d4bb536fSBarry Smith } else { 1396a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option"); 1397d4bb536fSBarry Smith } 13985a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1399c1f60f51SBarry Smith } 1400272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);CHKERRQ(ierr); 1401c1f60f51SBarry Smith /* 1402dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1403c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1404c1f60f51SBarry Smith */ 140531615d04SLois Curfman McInnes if (flg) { 1406272ac6f2SLois Curfman McInnes Mat J; 14075a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1408272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1409272ac6f2SLois Curfman McInnes snes->mfshell = J; 141083e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 141183e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP);CHKERRQ(ierr); 141294a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1413d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 141483e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP);CHKERRQ(ierr); 141594a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1416d4bb536fSBarry Smith } else { 1417a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option"); 1418d4bb536fSBarry Smith } 14195a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1420272ac6f2SLois Curfman McInnes } 14219b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 1422a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1423a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 14245a655dc6SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1425a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 1426a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector"); 1427a8c6a408SBarry Smith } 1428a703fe33SLois Curfman McInnes 1429a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1430*454a90a3SBarry Smith if (snes->ksp_ewconv && !PetscTypeCompare(snes,SNES_EQ_TR)) { 1431a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1432a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 1433a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 14345a655dc6SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,(void *)snes);CHKERRQ(ierr); 1435a703fe33SLois Curfman McInnes } 14369b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1437a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1438a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1439a8c6a408SBarry Smith if (!snes->computeumfunction) { 1440a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first"); 1441a8c6a408SBarry Smith } 14425a655dc6SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian()"); 1443d4bb536fSBarry Smith } else { 1444a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class"); 1445d4bb536fSBarry Smith } 1446a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 144782bf6240SBarry Smith snes->setupcalled = 1; 14483a40ed3dSBarry Smith PetscFunctionReturn(0); 14499b94acceSBarry Smith } 14509b94acceSBarry Smith 14515615d1e5SSatish Balay #undef __FUNC__ 1452d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy" 14539b94acceSBarry Smith /*@C 14549b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 14559b94acceSBarry Smith with SNESCreate(). 14569b94acceSBarry Smith 1457c7afd0dbSLois Curfman McInnes Collective on SNES 1458c7afd0dbSLois Curfman McInnes 14599b94acceSBarry Smith Input Parameter: 14609b94acceSBarry Smith . snes - the SNES context 14619b94acceSBarry Smith 146236851e7fSLois Curfman McInnes Level: beginner 146336851e7fSLois Curfman McInnes 14649b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 14659b94acceSBarry Smith 146663a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 14679b94acceSBarry Smith @*/ 14689b94acceSBarry Smith int SNESDestroy(SNES snes) 14699b94acceSBarry Smith { 1470b8a78c4aSBarry Smith int i,ierr; 14713a40ed3dSBarry Smith 14723a40ed3dSBarry Smith PetscFunctionBegin; 147377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 14743a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1475d4bb536fSBarry Smith 1476be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 1477be0abb6dSBarry Smith ierr = PetscAMSDestroy(snes);CHKERRQ(ierr); 1478be0abb6dSBarry Smith 1479e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1480606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 1481522c5e43SBarry Smith if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);} 14829b94acceSBarry Smith ierr = SLESDestroy(snes->sles);CHKERRQ(ierr); 1483522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1484b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++ ) { 1485b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1486b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1487b8a78c4aSBarry Smith } 1488b8a78c4aSBarry Smith } 14899b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 14900452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 14913a40ed3dSBarry Smith PetscFunctionReturn(0); 14929b94acceSBarry Smith } 14939b94acceSBarry Smith 14949b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 14959b94acceSBarry Smith 14965615d1e5SSatish Balay #undef __FUNC__ 14975615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 14989b94acceSBarry Smith /*@ 1499d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 15009b94acceSBarry Smith 1501c7afd0dbSLois Curfman McInnes Collective on SNES 1502c7afd0dbSLois Curfman McInnes 15039b94acceSBarry Smith Input Parameters: 1504c7afd0dbSLois Curfman McInnes + snes - the SNES context 150533174efeSLois Curfman McInnes . atol - absolute convergence tolerance 150633174efeSLois Curfman McInnes . rtol - relative convergence tolerance 150733174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 150833174efeSLois Curfman McInnes of the change in the solution between steps 150933174efeSLois Curfman McInnes . maxit - maximum number of iterations 1510c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1511fee21e36SBarry Smith 151233174efeSLois Curfman McInnes Options Database Keys: 1513c7afd0dbSLois Curfman McInnes + -snes_atol <atol> - Sets atol 1514c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1515c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1516c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1517c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 15189b94acceSBarry Smith 1519d7a720efSLois Curfman McInnes Notes: 15209b94acceSBarry Smith The default maximum number of iterations is 50. 15219b94acceSBarry Smith The default maximum number of function evaluations is 1000. 15229b94acceSBarry Smith 152336851e7fSLois Curfman McInnes Level: intermediate 152436851e7fSLois Curfman McInnes 152533174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 15269b94acceSBarry Smith 1527d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 15289b94acceSBarry Smith @*/ 1529d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 15309b94acceSBarry Smith { 15313a40ed3dSBarry Smith PetscFunctionBegin; 153277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1533d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1534d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1535d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1536d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1537d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 15383a40ed3dSBarry Smith PetscFunctionReturn(0); 15399b94acceSBarry Smith } 15409b94acceSBarry Smith 15415615d1e5SSatish Balay #undef __FUNC__ 15425615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 15439b94acceSBarry Smith /*@ 154433174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 154533174efeSLois Curfman McInnes 1546c7afd0dbSLois Curfman McInnes Not Collective 1547c7afd0dbSLois Curfman McInnes 154833174efeSLois Curfman McInnes Input Parameters: 1549c7afd0dbSLois Curfman McInnes + snes - the SNES context 155033174efeSLois Curfman McInnes . atol - absolute convergence tolerance 155133174efeSLois Curfman McInnes . rtol - relative convergence tolerance 155233174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 155333174efeSLois Curfman McInnes of the change in the solution between steps 155433174efeSLois Curfman McInnes . maxit - maximum number of iterations 1555c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1556fee21e36SBarry Smith 155733174efeSLois Curfman McInnes Notes: 155833174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 155933174efeSLois Curfman McInnes 156036851e7fSLois Curfman McInnes Level: intermediate 156136851e7fSLois Curfman McInnes 156233174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 156333174efeSLois Curfman McInnes 156433174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 156533174efeSLois Curfman McInnes @*/ 156633174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 156733174efeSLois Curfman McInnes { 15683a40ed3dSBarry Smith PetscFunctionBegin; 156933174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 157033174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 157133174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 157233174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 157333174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 157433174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 15753a40ed3dSBarry Smith PetscFunctionReturn(0); 157633174efeSLois Curfman McInnes } 157733174efeSLois Curfman McInnes 15785615d1e5SSatish Balay #undef __FUNC__ 15795615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 158033174efeSLois Curfman McInnes /*@ 15819b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 15829b94acceSBarry Smith 1583fee21e36SBarry Smith Collective on SNES 1584fee21e36SBarry Smith 1585c7afd0dbSLois Curfman McInnes Input Parameters: 1586c7afd0dbSLois Curfman McInnes + snes - the SNES context 1587c7afd0dbSLois Curfman McInnes - tol - tolerance 1588c7afd0dbSLois Curfman McInnes 15899b94acceSBarry Smith Options Database Key: 1590c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 15919b94acceSBarry Smith 159236851e7fSLois Curfman McInnes Level: intermediate 159336851e7fSLois Curfman McInnes 15949b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 15959b94acceSBarry Smith 1596d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 15979b94acceSBarry Smith @*/ 15989b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 15999b94acceSBarry Smith { 16003a40ed3dSBarry Smith PetscFunctionBegin; 160177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16029b94acceSBarry Smith snes->deltatol = tol; 16033a40ed3dSBarry Smith PetscFunctionReturn(0); 16049b94acceSBarry Smith } 16059b94acceSBarry Smith 16065615d1e5SSatish Balay #undef __FUNC__ 16075615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 16089b94acceSBarry Smith /*@ 160977c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 16109b94acceSBarry Smith for unconstrained minimization solvers. 16119b94acceSBarry Smith 1612fee21e36SBarry Smith Collective on SNES 1613fee21e36SBarry Smith 1614c7afd0dbSLois Curfman McInnes Input Parameters: 1615c7afd0dbSLois Curfman McInnes + snes - the SNES context 1616c7afd0dbSLois Curfman McInnes - ftol - minimum function tolerance 1617c7afd0dbSLois Curfman McInnes 16189b94acceSBarry Smith Options Database Key: 1619c7afd0dbSLois Curfman McInnes . -snes_fmin <ftol> - Sets ftol 16209b94acceSBarry Smith 16219b94acceSBarry Smith Note: 162277c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 16239b94acceSBarry Smith methods only. 16249b94acceSBarry Smith 162536851e7fSLois Curfman McInnes Level: intermediate 162636851e7fSLois Curfman McInnes 16279b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 16289b94acceSBarry Smith 1629d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 16309b94acceSBarry Smith @*/ 163177c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 16329b94acceSBarry Smith { 16333a40ed3dSBarry Smith PetscFunctionBegin; 163477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16359b94acceSBarry Smith snes->fmin = ftol; 16363a40ed3dSBarry Smith PetscFunctionReturn(0); 16379b94acceSBarry Smith } 1638df9fa365SBarry Smith /* 1639df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1640df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1641df9fa365SBarry Smith macros instead of functions 1642df9fa365SBarry Smith */ 1643ce1608b8SBarry Smith #undef __FUNC__ 1644ce1608b8SBarry Smith #define __FUNC__ "SNESLGMonitor" 1645ce1608b8SBarry Smith int SNESLGMonitor(SNES snes,int it,double norm,void *ctx) 1646ce1608b8SBarry Smith { 1647ce1608b8SBarry Smith int ierr; 1648ce1608b8SBarry Smith 1649ce1608b8SBarry Smith PetscFunctionBegin; 1650184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1651ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1652ce1608b8SBarry Smith PetscFunctionReturn(0); 1653ce1608b8SBarry Smith } 1654ce1608b8SBarry Smith 1655df9fa365SBarry Smith #undef __FUNC__ 1656df9fa365SBarry Smith #define __FUNC__ "SNESLGMonitorCreate" 1657df9fa365SBarry Smith int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n, DrawLG *draw) 1658df9fa365SBarry Smith { 1659df9fa365SBarry Smith int ierr; 1660df9fa365SBarry Smith 1661df9fa365SBarry Smith PetscFunctionBegin; 1662df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1663df9fa365SBarry Smith PetscFunctionReturn(0); 1664df9fa365SBarry Smith } 1665df9fa365SBarry Smith 1666df9fa365SBarry Smith #undef __FUNC__ 1667df9fa365SBarry Smith #define __FUNC__ "SNESLGMonitorDestroy" 1668df9fa365SBarry Smith int SNESLGMonitorDestroy(DrawLG draw) 1669df9fa365SBarry Smith { 1670df9fa365SBarry Smith int ierr; 1671df9fa365SBarry Smith 1672df9fa365SBarry Smith PetscFunctionBegin; 1673df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1674df9fa365SBarry Smith PetscFunctionReturn(0); 1675df9fa365SBarry Smith } 1676df9fa365SBarry Smith 16779b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 16789b94acceSBarry Smith 16795615d1e5SSatish Balay #undef __FUNC__ 1680d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor" 16819b94acceSBarry Smith /*@C 16825cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 16839b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 16849b94acceSBarry Smith progress. 16859b94acceSBarry Smith 1686fee21e36SBarry Smith Collective on SNES 1687fee21e36SBarry Smith 1688c7afd0dbSLois Curfman McInnes Input Parameters: 1689c7afd0dbSLois Curfman McInnes + snes - the SNES context 1690c7afd0dbSLois Curfman McInnes . func - monitoring routine 1691b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1692c7afd0dbSLois Curfman McInnes monitor routine (may be PETSC_NULL) 1693b8a78c4aSBarry Smith - monitordestroy - options routine that frees monitor context 16949b94acceSBarry Smith 1695c7afd0dbSLois Curfman McInnes Calling sequence of func: 169640a0c1c6SLois Curfman McInnes $ int func(SNES snes,int its, double norm,void *mctx) 1697c7afd0dbSLois Curfman McInnes 1698c7afd0dbSLois Curfman McInnes + snes - the SNES context 1699c7afd0dbSLois Curfman McInnes . its - iteration number 1700c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 1701c7afd0dbSLois Curfman McInnes for SNES_NONLINEAR_EQUATIONS methods 170240a0c1c6SLois Curfman McInnes . norm - 2-norm gradient value (may be estimated) 1703c7afd0dbSLois Curfman McInnes for SNES_UNCONSTRAINED_MINIMIZATION methods 170440a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 17059b94acceSBarry Smith 17069665c990SLois Curfman McInnes Options Database Keys: 1707c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1708c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1709c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1710c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1711c7afd0dbSLois Curfman McInnes been hardwired into a code by 1712c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1713c7afd0dbSLois Curfman McInnes does not cancel those set via 1714c7afd0dbSLois Curfman McInnes the options database. 17159665c990SLois Curfman McInnes 1716639f9d9dSBarry Smith Notes: 17176bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 17186bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 17196bc08f3fSLois Curfman McInnes order in which they were set. 1720639f9d9dSBarry Smith 172136851e7fSLois Curfman McInnes Level: intermediate 172236851e7fSLois Curfman McInnes 17239b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 17249b94acceSBarry Smith 17255cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 17269b94acceSBarry Smith @*/ 1727b8a78c4aSBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx,int (*monitordestroy)(void *)) 17289b94acceSBarry Smith { 17293a40ed3dSBarry Smith PetscFunctionBegin; 1730184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1731639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1732a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set"); 1733639f9d9dSBarry Smith } 1734639f9d9dSBarry Smith 1735639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1736b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1737639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 17383a40ed3dSBarry Smith PetscFunctionReturn(0); 17399b94acceSBarry Smith } 17409b94acceSBarry Smith 17415615d1e5SSatish Balay #undef __FUNC__ 17425cd90555SBarry Smith #define __FUNC__ "SNESClearMonitor" 17435cd90555SBarry Smith /*@C 17445cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 17455cd90555SBarry Smith 1746c7afd0dbSLois Curfman McInnes Collective on SNES 1747c7afd0dbSLois Curfman McInnes 17485cd90555SBarry Smith Input Parameters: 17495cd90555SBarry Smith . snes - the SNES context 17505cd90555SBarry Smith 17515cd90555SBarry Smith Options Database: 1752c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1753c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1754c7afd0dbSLois Curfman McInnes set via the options database 17555cd90555SBarry Smith 17565cd90555SBarry Smith Notes: 17575cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 17585cd90555SBarry Smith 175936851e7fSLois Curfman McInnes Level: intermediate 176036851e7fSLois Curfman McInnes 17615cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 17625cd90555SBarry Smith 17635cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 17645cd90555SBarry Smith @*/ 17655cd90555SBarry Smith int SNESClearMonitor( SNES snes ) 17665cd90555SBarry Smith { 17675cd90555SBarry Smith PetscFunctionBegin; 1768184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17695cd90555SBarry Smith snes->numbermonitors = 0; 17705cd90555SBarry Smith PetscFunctionReturn(0); 17715cd90555SBarry Smith } 17725cd90555SBarry Smith 17735cd90555SBarry Smith #undef __FUNC__ 1774d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest" 17759b94acceSBarry Smith /*@C 17769b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 17779b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 17789b94acceSBarry Smith 1779fee21e36SBarry Smith Collective on SNES 1780fee21e36SBarry Smith 1781c7afd0dbSLois Curfman McInnes Input Parameters: 1782c7afd0dbSLois Curfman McInnes + snes - the SNES context 1783c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1784c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1785c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 17869b94acceSBarry Smith 1787c7afd0dbSLois Curfman McInnes Calling sequence of func: 1788184914b5SBarry Smith $ int func (SNES snes,double xnorm,double gnorm,double f,SNESConvergedReason *reason,void *cctx) 1789c7afd0dbSLois Curfman McInnes 1790c7afd0dbSLois Curfman McInnes + snes - the SNES context 1791c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1792184914b5SBarry Smith . reason - reason for convergence/divergence 1793c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 1794c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods) 1795c7afd0dbSLois Curfman McInnes . f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods) 1796c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods) 1797c7afd0dbSLois Curfman McInnes - f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods) 17989b94acceSBarry Smith 179936851e7fSLois Curfman McInnes Level: advanced 180036851e7fSLois Curfman McInnes 18019b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 18029b94acceSBarry Smith 180340191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 180440191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 18059b94acceSBarry Smith @*/ 1806184914b5SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,SNESConvergedReason*,void*),void *cctx) 18079b94acceSBarry Smith { 18083a40ed3dSBarry Smith PetscFunctionBegin; 1809184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18109b94acceSBarry Smith (snes)->converged = func; 18119b94acceSBarry Smith (snes)->cnvP = cctx; 18123a40ed3dSBarry Smith PetscFunctionReturn(0); 18139b94acceSBarry Smith } 18149b94acceSBarry Smith 18155615d1e5SSatish Balay #undef __FUNC__ 1816184914b5SBarry Smith #define __FUNC__ "SNESGetConvergedReason" 1817184914b5SBarry Smith /*@C 1818184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1819184914b5SBarry Smith 1820184914b5SBarry Smith Not Collective 1821184914b5SBarry Smith 1822184914b5SBarry Smith Input Parameter: 1823184914b5SBarry Smith . snes - the SNES context 1824184914b5SBarry Smith 1825184914b5SBarry Smith Output Parameter: 1826184914b5SBarry Smith . reason - negative value indicates diverged, positive value converged, see snes.h or the 1827184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1828184914b5SBarry Smith 1829184914b5SBarry Smith Level: intermediate 1830184914b5SBarry Smith 1831184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1832184914b5SBarry Smith 1833184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1834184914b5SBarry Smith 1835184914b5SBarry Smith .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1836184914b5SBarry Smith SNESConverged_UM_LS(), SNESConverged_UM_TR() 1837184914b5SBarry Smith @*/ 1838184914b5SBarry Smith int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1839184914b5SBarry Smith { 1840184914b5SBarry Smith PetscFunctionBegin; 1841184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1842184914b5SBarry Smith *reason = snes->reason; 1843184914b5SBarry Smith PetscFunctionReturn(0); 1844184914b5SBarry Smith } 1845184914b5SBarry Smith 1846184914b5SBarry Smith #undef __FUNC__ 1847d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" 1848c9005455SLois Curfman McInnes /*@ 1849c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1850c9005455SLois Curfman McInnes 1851fee21e36SBarry Smith Collective on SNES 1852fee21e36SBarry Smith 1853c7afd0dbSLois Curfman McInnes Input Parameters: 1854c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1855c7afd0dbSLois Curfman McInnes . a - array to hold history 1856758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1857758f92a0SBarry Smith negative if not converged) for each solve. 1858758f92a0SBarry Smith . na - size of a and its 1859758f92a0SBarry Smith - reset - PETSC_TRUTH indicates each new nonlinear solve resets the history counter to zero, 1860758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1861c7afd0dbSLois Curfman McInnes 1862c9005455SLois Curfman McInnes Notes: 1863c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1864c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1865c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1866c9005455SLois Curfman McInnes at each step. 1867c9005455SLois Curfman McInnes 1868c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1869c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1870c9005455SLois Curfman McInnes during the section of code that is being timed. 1871c9005455SLois Curfman McInnes 187236851e7fSLois Curfman McInnes Level: intermediate 187336851e7fSLois Curfman McInnes 1874c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1875758f92a0SBarry Smith 187608405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1877758f92a0SBarry Smith 1878c9005455SLois Curfman McInnes @*/ 1879758f92a0SBarry Smith int SNESSetConvergenceHistory(SNES snes, double *a, int *its,int na,PetscTruth reset) 1880c9005455SLois Curfman McInnes { 18813a40ed3dSBarry Smith PetscFunctionBegin; 1882c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1883c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1884c9005455SLois Curfman McInnes snes->conv_hist = a; 1885758f92a0SBarry Smith snes->conv_hist_its = its; 1886758f92a0SBarry Smith snes->conv_hist_max = na; 1887758f92a0SBarry Smith snes->conv_hist_reset = reset; 1888758f92a0SBarry Smith PetscFunctionReturn(0); 1889758f92a0SBarry Smith } 1890758f92a0SBarry Smith 1891758f92a0SBarry Smith #undef __FUNC__ 1892758f92a0SBarry Smith #define __FUNC__ "SNESGetConvergenceHistory" 18930c4c9dddSBarry Smith /*@C 1894758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1895758f92a0SBarry Smith 1896758f92a0SBarry Smith Collective on SNES 1897758f92a0SBarry Smith 1898758f92a0SBarry Smith Input Parameter: 1899758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1900758f92a0SBarry Smith 1901758f92a0SBarry Smith Output Parameters: 1902758f92a0SBarry Smith . a - array to hold history 1903758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1904758f92a0SBarry Smith negative if not converged) for each solve. 1905758f92a0SBarry Smith - na - size of a and its 1906758f92a0SBarry Smith 1907758f92a0SBarry Smith Notes: 1908758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1909758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1910758f92a0SBarry Smith 1911758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1912758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1913758f92a0SBarry Smith during the section of code that is being timed. 1914758f92a0SBarry Smith 1915758f92a0SBarry Smith Level: intermediate 1916758f92a0SBarry Smith 1917758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1918758f92a0SBarry Smith 1919758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1920758f92a0SBarry Smith 1921758f92a0SBarry Smith @*/ 1922758f92a0SBarry Smith int SNESGetConvergenceHistory(SNES snes, double **a, int **its,int *na) 1923758f92a0SBarry Smith { 1924758f92a0SBarry Smith PetscFunctionBegin; 1925758f92a0SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1926758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1927758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1928758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 19293a40ed3dSBarry Smith PetscFunctionReturn(0); 1930c9005455SLois Curfman McInnes } 1931c9005455SLois Curfman McInnes 1932c9005455SLois Curfman McInnes #undef __FUNC__ 19335615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 19349b94acceSBarry Smith /* 19359b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 19369b94acceSBarry Smith positive parameter delta. 19379b94acceSBarry Smith 19389b94acceSBarry Smith Input Parameters: 1939c7afd0dbSLois Curfman McInnes + snes - the SNES context 19409b94acceSBarry Smith . y - approximate solution of linear system 19419b94acceSBarry Smith . fnorm - 2-norm of current function 1942c7afd0dbSLois Curfman McInnes - delta - trust region size 19439b94acceSBarry Smith 19449b94acceSBarry Smith Output Parameters: 1945c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 19469b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 19479b94acceSBarry Smith region, and exceeds zero otherwise. 1948c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 19499b94acceSBarry Smith 19509b94acceSBarry Smith Note: 195140191667SLois Curfman McInnes For non-trust region methods such as SNES_EQ_LS, the parameter delta 19529b94acceSBarry Smith is set to be the maximum allowable step size. 19539b94acceSBarry Smith 19549b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 19559b94acceSBarry Smith */ 19569b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 19579b94acceSBarry Smith double *gpnorm,double *ynorm) 19589b94acceSBarry Smith { 19599b94acceSBarry Smith double norm; 19609b94acceSBarry Smith Scalar cnorm; 19613a40ed3dSBarry Smith int ierr; 19623a40ed3dSBarry Smith 19633a40ed3dSBarry Smith PetscFunctionBegin; 1964184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1965184914b5SBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE); 1966184914b5SBarry Smith PetscCheckSameComm(snes,y); 1967184914b5SBarry Smith 19683a40ed3dSBarry Smith ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr); 19699b94acceSBarry Smith if (norm > *delta) { 19709b94acceSBarry Smith norm = *delta/norm; 19719b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 19729b94acceSBarry Smith cnorm = norm; 19739b94acceSBarry Smith VecScale( &cnorm, y ); 19749b94acceSBarry Smith *ynorm = *delta; 19759b94acceSBarry Smith } else { 19769b94acceSBarry Smith *gpnorm = 0.0; 19779b94acceSBarry Smith *ynorm = norm; 19789b94acceSBarry Smith } 19793a40ed3dSBarry Smith PetscFunctionReturn(0); 19809b94acceSBarry Smith } 19819b94acceSBarry Smith 19825615d1e5SSatish Balay #undef __FUNC__ 19835615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 19849b94acceSBarry Smith /*@ 19859b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 198663a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 19879b94acceSBarry Smith 1988c7afd0dbSLois Curfman McInnes Collective on SNES 1989c7afd0dbSLois Curfman McInnes 1990b2002411SLois Curfman McInnes Input Parameters: 1991c7afd0dbSLois Curfman McInnes + snes - the SNES context 1992c7afd0dbSLois Curfman McInnes - x - the solution vector 19939b94acceSBarry Smith 19949b94acceSBarry Smith Output Parameter: 1995b2002411SLois Curfman McInnes . its - number of iterations until termination 19969b94acceSBarry Smith 1997b2002411SLois Curfman McInnes Notes: 19988ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 19998ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 20008ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 20018ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 20028ddd3da0SLois Curfman McInnes 200336851e7fSLois Curfman McInnes Level: beginner 200436851e7fSLois Curfman McInnes 20059b94acceSBarry Smith .keywords: SNES, nonlinear, solve 20069b94acceSBarry Smith 200763a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 20089b94acceSBarry Smith @*/ 20098ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 20109b94acceSBarry Smith { 20113c7409f5SSatish Balay int ierr, flg; 2012052efed2SBarry Smith 20133a40ed3dSBarry Smith PetscFunctionBegin; 201477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2015184914b5SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 2016184914b5SBarry Smith PetscCheckSameComm(snes,x); 201774679c65SBarry Smith PetscValidIntPointer(its); 201882bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 2019c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 2020758f92a0SBarry Smith if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 20219b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 2022c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 20239b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr); 20249b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 20253c7409f5SSatish Balay ierr = OptionsHasName(PETSC_NULL,"-snes_view", &flg);CHKERRQ(ierr); 20266d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD);CHKERRQ(ierr); } 20273a40ed3dSBarry Smith PetscFunctionReturn(0); 20289b94acceSBarry Smith } 20299b94acceSBarry Smith 20309b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 20319b94acceSBarry Smith 20325615d1e5SSatish Balay #undef __FUNC__ 20335615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 203482bf6240SBarry Smith /*@C 20354b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 20369b94acceSBarry Smith 2037fee21e36SBarry Smith Collective on SNES 2038fee21e36SBarry Smith 2039c7afd0dbSLois Curfman McInnes Input Parameters: 2040c7afd0dbSLois Curfman McInnes + snes - the SNES context 2041*454a90a3SBarry Smith - type - a known method 2042c7afd0dbSLois Curfman McInnes 2043c7afd0dbSLois Curfman McInnes Options Database Key: 2044*454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 2045c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 2046ae12b187SLois Curfman McInnes 20479b94acceSBarry Smith Notes: 20489b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 204936851e7fSLois Curfman McInnes + SNES_EQ_LS - Newton's method with line search 2050c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 2051c7afd0dbSLois Curfman McInnes . SNES_EQ_TR - Newton's method with trust region 2052c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 2053c7afd0dbSLois Curfman McInnes . SNES_UM_TR - Newton's method with trust region 2054c7afd0dbSLois Curfman McInnes (unconstrained minimization) 205536851e7fSLois Curfman McInnes - SNES_UM_LS - Newton's method with line search 2056c7afd0dbSLois Curfman McInnes (unconstrained minimization) 20579b94acceSBarry Smith 2058ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 2059ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 2060ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 2061ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 2062ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 2063ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 2064ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 2065ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 2066ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 206736851e7fSLois Curfman McInnes appropriate method. In other words, this routine is not for beginners. 206836851e7fSLois Curfman McInnes 206936851e7fSLois Curfman McInnes Level: intermediate 2070a703fe33SLois Curfman McInnes 2071*454a90a3SBarry Smith .keywords: SNES, set, type 20729b94acceSBarry Smith @*/ 2073*454a90a3SBarry Smith int SNESSetType(SNES snes,SNESType type) 20749b94acceSBarry Smith { 207584cb2905SBarry Smith int ierr; 20769b94acceSBarry Smith int (*r)(SNES); 20773a40ed3dSBarry Smith 20783a40ed3dSBarry Smith PetscFunctionBegin; 207977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 208082bf6240SBarry Smith 2081*454a90a3SBarry Smith if (PetscTypeCompare(snes,type)) PetscFunctionReturn(0); 208282bf6240SBarry Smith 208382bf6240SBarry Smith if (snes->setupcalled) { 2084e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 208582bf6240SBarry Smith snes->data = 0; 208682bf6240SBarry Smith } 20877d1a2b2bSBarry Smith 20889b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 208982bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 209082bf6240SBarry Smith 2091*454a90a3SBarry Smith ierr = FListFind(snes->comm, SNESList, type,(int (**)(void *)) &r );CHKERRQ(ierr); 209282bf6240SBarry Smith 2093*454a90a3SBarry Smith if (!r) SETERRQ1(1,1,"Unable to find requested SNES type %s",type); 20941548b14aSBarry Smith 2095606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 209682bf6240SBarry Smith snes->data = 0; 20973a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 209882bf6240SBarry Smith 2099*454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 210082bf6240SBarry Smith snes->set_method_called = 1; 210182bf6240SBarry Smith 21023a40ed3dSBarry Smith PetscFunctionReturn(0); 21039b94acceSBarry Smith } 21049b94acceSBarry Smith 2105a847f771SSatish Balay 21069b94acceSBarry Smith /* --------------------------------------------------------------------- */ 21075615d1e5SSatish Balay #undef __FUNC__ 2108d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy" 21099b94acceSBarry Smith /*@C 21109b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 21119b94acceSBarry Smith registered by SNESRegister(). 21129b94acceSBarry Smith 2113fee21e36SBarry Smith Not Collective 2114fee21e36SBarry Smith 211536851e7fSLois Curfman McInnes Level: advanced 211636851e7fSLois Curfman McInnes 21179b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 21189b94acceSBarry Smith 21199b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 21209b94acceSBarry Smith @*/ 2121cf256101SBarry Smith int SNESRegisterDestroy(void) 21229b94acceSBarry Smith { 212382bf6240SBarry Smith int ierr; 212482bf6240SBarry Smith 21253a40ed3dSBarry Smith PetscFunctionBegin; 212682bf6240SBarry Smith if (SNESList) { 2127488ecbafSBarry Smith ierr = FListDestroy( SNESList );CHKERRQ(ierr); 212882bf6240SBarry Smith SNESList = 0; 21299b94acceSBarry Smith } 213084cb2905SBarry Smith SNESRegisterAllCalled = 0; 21313a40ed3dSBarry Smith PetscFunctionReturn(0); 21329b94acceSBarry Smith } 21339b94acceSBarry Smith 21345615d1e5SSatish Balay #undef __FUNC__ 2135d4bb536fSBarry Smith #define __FUNC__ "SNESGetType" 21369b94acceSBarry Smith /*@C 21379a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 21389b94acceSBarry Smith 2139c7afd0dbSLois Curfman McInnes Not Collective 2140c7afd0dbSLois Curfman McInnes 21419b94acceSBarry Smith Input Parameter: 21424b0e389bSBarry Smith . snes - nonlinear solver context 21439b94acceSBarry Smith 21449b94acceSBarry Smith Output Parameter: 2145*454a90a3SBarry Smith . type - SNES method (a charactor string) 21469b94acceSBarry Smith 214736851e7fSLois Curfman McInnes Level: intermediate 214836851e7fSLois Curfman McInnes 2149*454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 21509b94acceSBarry Smith @*/ 2151*454a90a3SBarry Smith int SNESGetType(SNES snes, SNESType *type) 21529b94acceSBarry Smith { 21533a40ed3dSBarry Smith PetscFunctionBegin; 2154184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2155*454a90a3SBarry Smith *type = snes->type_name; 21563a40ed3dSBarry Smith PetscFunctionReturn(0); 21579b94acceSBarry Smith } 21589b94acceSBarry Smith 21595615d1e5SSatish Balay #undef __FUNC__ 2160d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution" 21619b94acceSBarry Smith /*@C 21629b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 21639b94acceSBarry Smith stored. 21649b94acceSBarry Smith 2165c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2166c7afd0dbSLois Curfman McInnes 21679b94acceSBarry Smith Input Parameter: 21689b94acceSBarry Smith . snes - the SNES context 21699b94acceSBarry Smith 21709b94acceSBarry Smith Output Parameter: 21719b94acceSBarry Smith . x - the solution 21729b94acceSBarry Smith 217336851e7fSLois Curfman McInnes Level: advanced 217436851e7fSLois Curfman McInnes 21759b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 21769b94acceSBarry Smith 2177a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 21789b94acceSBarry Smith @*/ 21799b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 21809b94acceSBarry Smith { 21813a40ed3dSBarry Smith PetscFunctionBegin; 218277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 21839b94acceSBarry Smith *x = snes->vec_sol_always; 21843a40ed3dSBarry Smith PetscFunctionReturn(0); 21859b94acceSBarry Smith } 21869b94acceSBarry Smith 21875615d1e5SSatish Balay #undef __FUNC__ 2188d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" 21899b94acceSBarry Smith /*@C 21909b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 21919b94acceSBarry Smith stored. 21929b94acceSBarry Smith 2193c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2194c7afd0dbSLois Curfman McInnes 21959b94acceSBarry Smith Input Parameter: 21969b94acceSBarry Smith . snes - the SNES context 21979b94acceSBarry Smith 21989b94acceSBarry Smith Output Parameter: 21999b94acceSBarry Smith . x - the solution update 22009b94acceSBarry Smith 220136851e7fSLois Curfman McInnes Level: advanced 220236851e7fSLois Curfman McInnes 22039b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 22049b94acceSBarry Smith 22059b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 22069b94acceSBarry Smith @*/ 22079b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 22089b94acceSBarry Smith { 22093a40ed3dSBarry Smith PetscFunctionBegin; 221077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 22119b94acceSBarry Smith *x = snes->vec_sol_update_always; 22123a40ed3dSBarry Smith PetscFunctionReturn(0); 22139b94acceSBarry Smith } 22149b94acceSBarry Smith 22155615d1e5SSatish Balay #undef __FUNC__ 2216d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction" 22179b94acceSBarry Smith /*@C 22183638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 22199b94acceSBarry Smith 2220c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2221c7afd0dbSLois Curfman McInnes 22229b94acceSBarry Smith Input Parameter: 22239b94acceSBarry Smith . snes - the SNES context 22249b94acceSBarry Smith 22259b94acceSBarry Smith Output Parameter: 22267bf4e008SBarry Smith + r - the function (or PETSC_NULL) 22277bf4e008SBarry Smith - ctx - the function context (or PETSC_NULL) 22289b94acceSBarry Smith 22299b94acceSBarry Smith Notes: 22309b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 22319b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 22329b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 22339b94acceSBarry Smith 223436851e7fSLois Curfman McInnes Level: advanced 223536851e7fSLois Curfman McInnes 2236a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 22379b94acceSBarry Smith 223831615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 223931615d04SLois Curfman McInnes SNESGetGradient() 22407bf4e008SBarry Smith 22419b94acceSBarry Smith @*/ 22427bf4e008SBarry Smith int SNESGetFunction(SNES snes,Vec *r,void **ctx) 22439b94acceSBarry Smith { 22443a40ed3dSBarry Smith PetscFunctionBegin; 224577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2246a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 2247a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 2248a8c6a408SBarry Smith } 22497bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 22507bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 22513a40ed3dSBarry Smith PetscFunctionReturn(0); 22529b94acceSBarry Smith } 22539b94acceSBarry Smith 22545615d1e5SSatish Balay #undef __FUNC__ 2255d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient" 22569b94acceSBarry Smith /*@C 22573638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 22589b94acceSBarry Smith 2259c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2260c7afd0dbSLois Curfman McInnes 22619b94acceSBarry Smith Input Parameter: 22629b94acceSBarry Smith . snes - the SNES context 22639b94acceSBarry Smith 22649b94acceSBarry Smith Output Parameter: 22657bf4e008SBarry Smith + r - the gradient (or PETSC_NULL) 22667bf4e008SBarry Smith - ctx - the gradient context (or PETSC_NULL) 22679b94acceSBarry Smith 22689b94acceSBarry Smith Notes: 22699b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 22709b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 22719b94acceSBarry Smith SNESGetFunction(). 22729b94acceSBarry Smith 227336851e7fSLois Curfman McInnes Level: advanced 227436851e7fSLois Curfman McInnes 22759b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 22769b94acceSBarry Smith 22777bf4e008SBarry Smith .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction(), 22787bf4e008SBarry Smith SNESSetGradient(), SNESSetFunction() 22797bf4e008SBarry Smith 22809b94acceSBarry Smith @*/ 22817bf4e008SBarry Smith int SNESGetGradient(SNES snes,Vec *r,void **ctx) 22829b94acceSBarry Smith { 22833a40ed3dSBarry Smith PetscFunctionBegin; 228477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 22853a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 2286a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 22873a40ed3dSBarry Smith } 22887bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 22897bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 22903a40ed3dSBarry Smith PetscFunctionReturn(0); 22919b94acceSBarry Smith } 22929b94acceSBarry Smith 22935615d1e5SSatish Balay #undef __FUNC__ 2294d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" 22957bf4e008SBarry Smith /*@C 22969b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 22979b94acceSBarry Smith unconstrained minimization problems. 22989b94acceSBarry Smith 2299c7afd0dbSLois Curfman McInnes Not Collective 2300c7afd0dbSLois Curfman McInnes 23019b94acceSBarry Smith Input Parameter: 23029b94acceSBarry Smith . snes - the SNES context 23039b94acceSBarry Smith 23049b94acceSBarry Smith Output Parameter: 23057bf4e008SBarry Smith + r - the function (or PETSC_NULL) 23067bf4e008SBarry Smith - ctx - the function context (or PETSC_NULL) 23079b94acceSBarry Smith 23089b94acceSBarry Smith Notes: 23099b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 23109b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 23119b94acceSBarry Smith SNESGetFunction(). 23129b94acceSBarry Smith 231336851e7fSLois Curfman McInnes Level: advanced 231436851e7fSLois Curfman McInnes 23159b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 23169b94acceSBarry Smith 23177bf4e008SBarry Smith .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction(), SNESSetFunction() 23187bf4e008SBarry Smith 23199b94acceSBarry Smith @*/ 23207bf4e008SBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r,void **ctx) 23219b94acceSBarry Smith { 23223a40ed3dSBarry Smith PetscFunctionBegin; 232377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 232474679c65SBarry Smith PetscValidScalarPointer(r); 23253a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 2326a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 23273a40ed3dSBarry Smith } 23287bf4e008SBarry Smith if (r) *r = snes->fc; 23297bf4e008SBarry Smith if (ctx) *ctx = snes->umfunP; 23303a40ed3dSBarry Smith PetscFunctionReturn(0); 23319b94acceSBarry Smith } 23329b94acceSBarry Smith 23335615d1e5SSatish Balay #undef __FUNC__ 2334d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" 23353c7409f5SSatish Balay /*@C 23363c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 2337d850072dSLois Curfman McInnes SNES options in the database. 23383c7409f5SSatish Balay 2339fee21e36SBarry Smith Collective on SNES 2340fee21e36SBarry Smith 2341c7afd0dbSLois Curfman McInnes Input Parameter: 2342c7afd0dbSLois Curfman McInnes + snes - the SNES context 2343c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2344c7afd0dbSLois Curfman McInnes 2345d850072dSLois Curfman McInnes Notes: 2346a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2347c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2348d850072dSLois Curfman McInnes 234936851e7fSLois Curfman McInnes Level: advanced 235036851e7fSLois Curfman McInnes 23513c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 2352a86d99e1SLois Curfman McInnes 2353a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 23543c7409f5SSatish Balay @*/ 23553c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 23563c7409f5SSatish Balay { 23573c7409f5SSatish Balay int ierr; 23583c7409f5SSatish Balay 23593a40ed3dSBarry Smith PetscFunctionBegin; 236077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2361639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix);CHKERRQ(ierr); 23623c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 23633a40ed3dSBarry Smith PetscFunctionReturn(0); 23643c7409f5SSatish Balay } 23653c7409f5SSatish Balay 23665615d1e5SSatish Balay #undef __FUNC__ 2367d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" 23683c7409f5SSatish Balay /*@C 2369f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2370d850072dSLois Curfman McInnes SNES options in the database. 23713c7409f5SSatish Balay 2372fee21e36SBarry Smith Collective on SNES 2373fee21e36SBarry Smith 2374c7afd0dbSLois Curfman McInnes Input Parameters: 2375c7afd0dbSLois Curfman McInnes + snes - the SNES context 2376c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2377c7afd0dbSLois Curfman McInnes 2378d850072dSLois Curfman McInnes Notes: 2379a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2380c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2381d850072dSLois Curfman McInnes 238236851e7fSLois Curfman McInnes Level: advanced 238336851e7fSLois Curfman McInnes 23843c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2385a86d99e1SLois Curfman McInnes 2386a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 23873c7409f5SSatish Balay @*/ 23883c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 23893c7409f5SSatish Balay { 23903c7409f5SSatish Balay int ierr; 23913c7409f5SSatish Balay 23923a40ed3dSBarry Smith PetscFunctionBegin; 239377c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2394639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix);CHKERRQ(ierr); 23953c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 23963a40ed3dSBarry Smith PetscFunctionReturn(0); 23973c7409f5SSatish Balay } 23983c7409f5SSatish Balay 23995615d1e5SSatish Balay #undef __FUNC__ 2400d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" 24019ab63eb5SSatish Balay /*@C 24023c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 24033c7409f5SSatish Balay SNES options in the database. 24043c7409f5SSatish Balay 2405c7afd0dbSLois Curfman McInnes Not Collective 2406c7afd0dbSLois Curfman McInnes 24073c7409f5SSatish Balay Input Parameter: 24083c7409f5SSatish Balay . snes - the SNES context 24093c7409f5SSatish Balay 24103c7409f5SSatish Balay Output Parameter: 24113c7409f5SSatish Balay . prefix - pointer to the prefix string used 24123c7409f5SSatish Balay 24139ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 24149ab63eb5SSatish Balay sufficient length to hold the prefix. 24159ab63eb5SSatish Balay 241636851e7fSLois Curfman McInnes Level: advanced 241736851e7fSLois Curfman McInnes 24183c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2419a86d99e1SLois Curfman McInnes 2420a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 24213c7409f5SSatish Balay @*/ 24223c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 24233c7409f5SSatish Balay { 24243c7409f5SSatish Balay int ierr; 24253c7409f5SSatish Balay 24263a40ed3dSBarry Smith PetscFunctionBegin; 242777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2428639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix);CHKERRQ(ierr); 24293a40ed3dSBarry Smith PetscFunctionReturn(0); 24303c7409f5SSatish Balay } 24313c7409f5SSatish Balay 2432ca161407SBarry Smith #undef __FUNC__ 2433ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp" 2434ca161407SBarry Smith /*@ 2435ca161407SBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2436ca161407SBarry Smith 2437c7afd0dbSLois Curfman McInnes Collective on SNES 2438c7afd0dbSLois Curfman McInnes 2439ca161407SBarry Smith Input Parameter: 2440ca161407SBarry Smith . snes - the SNES context 2441ca161407SBarry Smith 2442ca161407SBarry Smith Options Database Keys: 2443c7afd0dbSLois Curfman McInnes + -help - Prints SNES options 2444c7afd0dbSLois Curfman McInnes - -h - Prints SNES options 2445ca161407SBarry Smith 244636851e7fSLois Curfman McInnes Level: beginner 244736851e7fSLois Curfman McInnes 2448ca161407SBarry Smith .keywords: SNES, nonlinear, help 2449ca161407SBarry Smith 2450ca161407SBarry Smith .seealso: SNESSetFromOptions() 2451ca161407SBarry Smith @*/ 2452ca161407SBarry Smith int SNESPrintHelp(SNES snes) 2453ca161407SBarry Smith { 2454ca161407SBarry Smith char p[64]; 2455ca161407SBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2456ca161407SBarry Smith int ierr; 2457ca161407SBarry Smith 2458ca161407SBarry Smith PetscFunctionBegin; 2459ca161407SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2460ca161407SBarry Smith 2461549d3d68SSatish Balay ierr = PetscStrcpy(p,"-");CHKERRQ(ierr); 2462ca161407SBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2463ca161407SBarry Smith 2464ca161407SBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2465ca161407SBarry Smith 2466ebb8b11fSBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 2467d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");CHKERRQ(ierr); 2468488ecbafSBarry Smith ierr = FListPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr); 2469d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);CHKERRQ(ierr); 2470d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);CHKERRQ(ierr); 2471d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);CHKERRQ(ierr); 2472d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);CHKERRQ(ierr); 2473d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);CHKERRQ(ierr); 2474d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);CHKERRQ(ierr); 2475d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);CHKERRQ(ierr); 2476d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");CHKERRQ(ierr); 2477d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);CHKERRQ(ierr); 2478d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 2479d132466eSBarry Smith residual norm at each iteration.\n",p);CHKERRQ(ierr); 2480d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 2481ca161407SBarry Smith residual norm for small residual norms. This is useful to conceal\n\ 2482d132466eSBarry Smith meaningless digits that may be different on different machines.\n",p);CHKERRQ(ierr); 2483d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);CHKERRQ(ierr); 2484d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_vecmonitor: plots solution at each iteration \n",p);CHKERRQ(ierr); 2485d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_vecmonitor_update: plots update to solution at each iteration \n",p);CHKERRQ(ierr); 2486ca161407SBarry Smith if (snes->type == SNES_NONLINEAR_EQUATIONS) { 2487d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2488d132466eSBarry Smith " Options for solving systems of nonlinear equations only:\n");CHKERRQ(ierr); 2489d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p);CHKERRQ(ierr); 2490d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p);CHKERRQ(ierr); 2491d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);CHKERRQ(ierr); 2492d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_mf_ksp_monitor - if using matrix-free multiply then prints h at each KSP iteration\n",p);CHKERRQ(ierr); 2493d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);CHKERRQ(ierr); 2494d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);CHKERRQ(ierr); 2495d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2496d132466eSBarry Smith " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);CHKERRQ(ierr); 2497d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2498d132466eSBarry Smith " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);CHKERRQ(ierr); 2499d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2500d132466eSBarry Smith " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);CHKERRQ(ierr); 2501d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2502d132466eSBarry Smith " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);CHKERRQ(ierr); 2503d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2504d132466eSBarry Smith " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);CHKERRQ(ierr); 2505d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2506d132466eSBarry Smith " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);CHKERRQ(ierr); 2507d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2508d132466eSBarry Smith " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);CHKERRQ(ierr); 2509ca161407SBarry Smith } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 2510d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," Options for solving unconstrained minimization problems only:\n");CHKERRQ(ierr); 2511d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);CHKERRQ(ierr); 2512d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p);CHKERRQ(ierr); 2513d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p);CHKERRQ(ierr); 2514ca161407SBarry Smith } 2515*454a90a3SBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <type> for help on ",p);CHKERRQ(ierr); 2516d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm,"a particular method\n");CHKERRQ(ierr); 2517ca161407SBarry Smith if (snes->printhelp) { 2518ca161407SBarry Smith ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr); 2519ca161407SBarry Smith } 2520ca161407SBarry Smith PetscFunctionReturn(0); 2521ca161407SBarry Smith } 25223c7409f5SSatish Balay 2523acb85ae6SSatish Balay /*MC 2524b2002411SLois Curfman McInnes SNESRegister - Adds a method to the nonlinear solver package. 25259b94acceSBarry Smith 2526b2002411SLois Curfman McInnes Synopsis: 2527b2002411SLois Curfman McInnes SNESRegister(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES)) 25289b94acceSBarry Smith 25298d76a1e5SLois Curfman McInnes Not collective 25308d76a1e5SLois Curfman McInnes 2531b2002411SLois Curfman McInnes Input Parameters: 2532c7afd0dbSLois Curfman McInnes + name_solver - name of a new user-defined solver 2533b2002411SLois Curfman McInnes . path - path (either absolute or relative) the library containing this solver 2534b2002411SLois Curfman McInnes . name_create - name of routine to create method context 2535c7afd0dbSLois Curfman McInnes - routine_create - routine to create method context 25369b94acceSBarry Smith 2537b2002411SLois Curfman McInnes Notes: 2538b2002411SLois Curfman McInnes SNESRegister() may be called multiple times to add several user-defined solvers. 25393a40ed3dSBarry Smith 2540b2002411SLois Curfman McInnes If dynamic libraries are used, then the fourth input argument (routine_create) 2541b2002411SLois Curfman McInnes is ignored. 2542b2002411SLois Curfman McInnes 2543b2002411SLois Curfman McInnes Sample usage: 254469225d78SLois Curfman McInnes .vb 2545b2002411SLois Curfman McInnes SNESRegister("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a, 2546b2002411SLois Curfman McInnes "MySolverCreate",MySolverCreate); 254769225d78SLois Curfman McInnes .ve 2548b2002411SLois Curfman McInnes 2549b2002411SLois Curfman McInnes Then, your solver can be chosen with the procedural interface via 2550b2002411SLois Curfman McInnes $ SNESSetType(snes,"my_solver") 2551b2002411SLois Curfman McInnes or at runtime via the option 2552b2002411SLois Curfman McInnes $ -snes_type my_solver 2553b2002411SLois Curfman McInnes 255436851e7fSLois Curfman McInnes Level: advanced 255536851e7fSLois Curfman McInnes 255685614651SBarry Smith $PETSC_ARCH, $PETSC_DIR, $PETSC_LDIR, and $BOPT occuring in pathname will be replaced with appropriate values. 2557dd438238SBarry Smith 2558b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register 2559b2002411SLois Curfman McInnes 2560b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy() 2561b2002411SLois Curfman McInnes M*/ 2562b2002411SLois Curfman McInnes 2563b2002411SLois Curfman McInnes #undef __FUNC__ 2564b2002411SLois Curfman McInnes #define __FUNC__ "SNESRegister_Private" 2565b2002411SLois Curfman McInnes int SNESRegister_Private(char *sname,char *path,char *name,int (*function)(SNES)) 2566b2002411SLois Curfman McInnes { 2567b2002411SLois Curfman McInnes char fullname[256]; 2568b2002411SLois Curfman McInnes int ierr; 2569b2002411SLois Curfman McInnes 2570b2002411SLois Curfman McInnes PetscFunctionBegin; 2571549d3d68SSatish Balay ierr = PetscStrcpy(fullname,path);CHKERRQ(ierr); 2572184914b5SBarry Smith ierr = PetscStrcat(fullname,":");CHKERRQ(ierr); 2573184914b5SBarry Smith ierr = PetscStrcat(fullname,name);CHKERRQ(ierr); 2574488ecbafSBarry Smith ierr = FListAdd_Private(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr); 2575b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2576b2002411SLois Curfman McInnes } 2577