1*c38d4ed2SBarry Smith /*$Id: snes.c,v 1.202 1999/11/05 14:47:05 bsmith Exp bsmith $*/ 29b94acceSBarry Smith 370f55243SBarry Smith #include "src/snes/snesimpl.h" /*I "snes.h" I*/ 49b94acceSBarry Smith 584cb2905SBarry Smith int SNESRegisterAllCalled = 0; 6488ecbafSBarry Smith FList SNESList = 0; 782bf6240SBarry Smith 85615d1e5SSatish Balay #undef __FUNC__ 9d4bb536fSBarry Smith #define __FUNC__ "SNESView" 109b94acceSBarry Smith /*@ 119b94acceSBarry Smith SNESView - Prints the SNES data structure. 129b94acceSBarry Smith 13fee21e36SBarry Smith Collective on SNES, unless Viewer is VIEWER_STDOUT_SELF 14fee21e36SBarry Smith 15c7afd0dbSLois Curfman McInnes Input Parameters: 16c7afd0dbSLois Curfman McInnes + SNES - the SNES context 17c7afd0dbSLois Curfman McInnes - viewer - visualization context 18c7afd0dbSLois Curfman McInnes 199b94acceSBarry Smith Options Database Key: 20c8a8ba5cSLois Curfman McInnes . -snes_view - Calls SNESView() at end of SNESSolve() 219b94acceSBarry Smith 229b94acceSBarry Smith Notes: 239b94acceSBarry Smith The available visualization contexts include 24c7afd0dbSLois Curfman McInnes + VIEWER_STDOUT_SELF - standard output (default) 25c7afd0dbSLois Curfman McInnes - VIEWER_STDOUT_WORLD - synchronized standard 26c8a8ba5cSLois Curfman McInnes output where only the first processor opens 27c8a8ba5cSLois Curfman McInnes the file. All other processors send their 28c8a8ba5cSLois Curfman McInnes data to the first processor to print. 299b94acceSBarry Smith 303e081fefSLois Curfman McInnes The user can open an alternative visualization context with 3177ed5343SBarry Smith ViewerASCIIOpen() - output to a specified file. 329b94acceSBarry Smith 3336851e7fSLois Curfman McInnes Level: beginner 3436851e7fSLois Curfman McInnes 359b94acceSBarry Smith .keywords: SNES, view 369b94acceSBarry Smith 3777ed5343SBarry Smith .seealso: ViewerASCIIOpen() 389b94acceSBarry Smith @*/ 399b94acceSBarry Smith int SNESView(SNES snes,Viewer viewer) 409b94acceSBarry Smith { 419b94acceSBarry Smith SNES_KSP_EW_ConvCtx *kctx; 429b94acceSBarry Smith int ierr; 439b94acceSBarry Smith SLES sles; 44454a90a3SBarry Smith char *type; 456831982aSBarry Smith PetscTruth isascii,isstring; 469b94acceSBarry Smith 473a40ed3dSBarry Smith PetscFunctionBegin; 4874679c65SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 490f5bd95cSBarry Smith if (!viewer) viewer = VIEWER_STDOUT_SELF; 500f5bd95cSBarry Smith PetscValidHeaderSpecific(viewer,VIEWER_COOKIE); 516831982aSBarry Smith PetscCheckSameComm(snes,viewer); 5274679c65SBarry Smith 536831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,ASCII_VIEWER,&isascii);CHKERRQ(ierr); 546831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)viewer,STRING_VIEWER,&isstring);CHKERRQ(ierr); 550f5bd95cSBarry Smith if (isascii) { 56d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 57454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 58454a90a3SBarry Smith if (type) { 59454a90a3SBarry Smith ierr = ViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 60184914b5SBarry Smith } else { 61454a90a3SBarry Smith ierr = ViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 62184914b5SBarry Smith } 630ef38995SBarry Smith if (snes->view) { 640ef38995SBarry Smith ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 650ef38995SBarry Smith ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 660ef38995SBarry Smith ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 670ef38995SBarry Smith } 68d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 69d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, truncation=%g, solution=%g\n", 70d132466eSBarry Smith snes->rtol, snes->atol, snes->trunctol, snes->xtol);CHKERRQ(ierr); 71d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr); 72d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr); 730ef38995SBarry Smith if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 74d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," min function tolerance=%g\n",snes->fmin);CHKERRQ(ierr); 750ef38995SBarry Smith } 769b94acceSBarry Smith if (snes->ksp_ewconv) { 779b94acceSBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 789b94acceSBarry Smith if (kctx) { 79d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr); 80d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 81d132466eSBarry Smith ierr = ViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 829b94acceSBarry Smith } 839b94acceSBarry Smith } 840f5bd95cSBarry Smith } else if (isstring) { 85454a90a3SBarry Smith ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 86454a90a3SBarry Smith ierr = ViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 8719bcc07fSBarry Smith } 8877ed5343SBarry Smith ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 890ef38995SBarry Smith ierr = ViewerASCIIPushTab(viewer);CHKERRQ(ierr); 909b94acceSBarry Smith ierr = SLESView(sles,viewer);CHKERRQ(ierr); 910ef38995SBarry Smith ierr = ViewerASCIIPopTab(viewer);CHKERRQ(ierr); 923a40ed3dSBarry Smith PetscFunctionReturn(0); 939b94acceSBarry Smith } 949b94acceSBarry Smith 95639f9d9dSBarry Smith /* 96639f9d9dSBarry Smith We retain a list of functions that also take SNES command 97639f9d9dSBarry Smith line options. These are called at the end SNESSetFromOptions() 98639f9d9dSBarry Smith */ 99639f9d9dSBarry Smith #define MAXSETFROMOPTIONS 5 100639f9d9dSBarry Smith static int numberofsetfromoptions; 101639f9d9dSBarry Smith static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 102639f9d9dSBarry Smith 1035615d1e5SSatish Balay #undef __FUNC__ 104d4bb536fSBarry Smith #define __FUNC__ "SNESAddOptionsChecker" 105639f9d9dSBarry Smith /*@ 106639f9d9dSBarry Smith SNESAddOptionsChecker - Adds an additional function to check for SNES options. 107639f9d9dSBarry Smith 108c7afd0dbSLois Curfman McInnes Not Collective 109c7afd0dbSLois Curfman McInnes 110639f9d9dSBarry Smith Input Parameter: 111639f9d9dSBarry Smith . snescheck - function that checks for options 112639f9d9dSBarry Smith 11336851e7fSLois Curfman McInnes Level: developer 11436851e7fSLois Curfman McInnes 115639f9d9dSBarry Smith .seealso: SNESSetFromOptions() 116639f9d9dSBarry Smith @*/ 117639f9d9dSBarry Smith int SNESAddOptionsChecker(int (*snescheck)(SNES) ) 118639f9d9dSBarry Smith { 1193a40ed3dSBarry Smith PetscFunctionBegin; 120639f9d9dSBarry Smith if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 121a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many options checkers, only 5 allowed"); 122639f9d9dSBarry Smith } 123639f9d9dSBarry Smith 124639f9d9dSBarry Smith othersetfromoptions[numberofsetfromoptions++] = snescheck; 1253a40ed3dSBarry Smith PetscFunctionReturn(0); 126639f9d9dSBarry Smith } 127639f9d9dSBarry Smith 1285615d1e5SSatish Balay #undef __FUNC__ 12915091d37SBarry Smith #define __FUNC__ "SNESSetTypeFromOptions" 13015091d37SBarry Smith /*@ 13115091d37SBarry Smith SNESSetTypeFromOptions - Sets the SNES solver type from the options database, 13215091d37SBarry Smith or sets a default if none is give. 13315091d37SBarry Smith 13415091d37SBarry Smith Collective on SNES 13515091d37SBarry Smith 13615091d37SBarry Smith Input Parameter: 13715091d37SBarry Smith . snes - the SNES context 13815091d37SBarry Smith 13915091d37SBarry Smith Options Database Keys: 1406831982aSBarry Smith . -snes_type <type> - ls, tr, umls, umtr, test 14115091d37SBarry Smith 14215091d37SBarry Smith Level: beginner 14315091d37SBarry Smith 14415091d37SBarry Smith .keywords: SNES, nonlinear, set, options, database 14515091d37SBarry Smith 14615091d37SBarry Smith .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetFromOptions() 14715091d37SBarry Smith @*/ 14815091d37SBarry Smith int SNESSetTypeFromOptions(SNES snes) 14915091d37SBarry Smith { 150454a90a3SBarry Smith char type[256]; 151f1af5d2fSBarry Smith int ierr; 152f1af5d2fSBarry Smith PetscTruth flg; 15315091d37SBarry Smith 15415091d37SBarry Smith PetscFunctionBegin; 15515091d37SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 15615091d37SBarry Smith if (snes->setupcalled) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call prior to SNESSetUp()"); 157454a90a3SBarry Smith ierr = OptionsGetString(snes->prefix,"-snes_type",type,256,&flg);CHKERRQ(ierr); 15815091d37SBarry Smith if (flg) { 159454a90a3SBarry Smith ierr = SNESSetType(snes,(SNESType) type);CHKERRQ(ierr); 16015091d37SBarry Smith } 16115091d37SBarry Smith /* 16215091d37SBarry Smith If SNES type has not yet been set, set it now 16315091d37SBarry Smith */ 16415091d37SBarry Smith if (!snes->type_name) { 16515091d37SBarry Smith if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1666831982aSBarry Smith ierr = SNESSetType(snes,SNESEQLS);CHKERRQ(ierr); 16715091d37SBarry Smith } else { 1686831982aSBarry Smith ierr = SNESSetType(snes,SNESUMTR);CHKERRQ(ierr); 16915091d37SBarry Smith } 17015091d37SBarry Smith } 17115091d37SBarry Smith PetscFunctionReturn(0); 17215091d37SBarry Smith } 17315091d37SBarry Smith 17415091d37SBarry Smith #undef __FUNC__ 1755615d1e5SSatish Balay #define __FUNC__ "SNESSetFromOptions" 1769b94acceSBarry Smith /*@ 177682d7d0cSBarry Smith SNESSetFromOptions - Sets various SNES and SLES parameters from user options. 1789b94acceSBarry Smith 179c7afd0dbSLois Curfman McInnes Collective on SNES 180c7afd0dbSLois Curfman McInnes 1819b94acceSBarry Smith Input Parameter: 1829b94acceSBarry Smith . snes - the SNES context 1839b94acceSBarry Smith 18436851e7fSLois Curfman McInnes Options Database Keys: 1856831982aSBarry Smith + -snes_type <type> - ls, tr, umls, umtr, test 18682738288SBarry Smith . -snes_stol - convergence tolerance in terms of the norm 18782738288SBarry Smith of the change in the solution between steps 188b39c3a46SLois Curfman McInnes . -snes_atol <atol> - absolute tolerance of residual norm 189b39c3a46SLois Curfman McInnes . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 190b39c3a46SLois Curfman McInnes . -snes_max_it <max_it> - maximum number of iterations 191b39c3a46SLois Curfman McInnes . -snes_max_funcs <max_funcs> - maximum number of function evaluations 192b39c3a46SLois Curfman McInnes . -snes_trtol <trtol> - trust region tolerance 19382738288SBarry Smith . -snes_no_convergence_test - skip convergence test in nonlinear or minimization 19482738288SBarry Smith solver; hence iterations will continue until max_it 1951fbbfb26SLois Curfman McInnes or some other criterion is reached. Saves expense 19682738288SBarry Smith of convergence test 19782738288SBarry Smith . -snes_monitor - prints residual norm at each iteration 198d132466eSBarry Smith . -snes_vecmonitor - plots solution at each iteration 199d132466eSBarry Smith . -snes_vecmonitor_update - plots update to solution at each iteration 20082738288SBarry Smith . -snes_xmonitor - plots residual norm at each iteration 201e24b481bSBarry Smith . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 20236851e7fSLois Curfman McInnes - -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 20382738288SBarry Smith 20482738288SBarry Smith Options Database for Eisenstat-Walker method: 20582738288SBarry Smith + -snes_ksp_eq_conv - use Eisenstat-Walker method for determining linear system convergence 20682738288SBarry Smith . -snes_ksp_eq_version ver - version of Eisenstat-Walker method 20736851e7fSLois Curfman McInnes . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 20836851e7fSLois Curfman McInnes . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 20936851e7fSLois Curfman McInnes . -snes_ksp_ew_gamma <gamma> - Sets gamma 21036851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha <alpha> - Sets alpha 21136851e7fSLois Curfman McInnes . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 21236851e7fSLois Curfman McInnes - -snes_ksp_ew_threshold <threshold> - Sets threshold 21382738288SBarry Smith 21411ca99fdSLois Curfman McInnes Notes: 21511ca99fdSLois Curfman McInnes To see all options, run your program with the -help option or consult 21611ca99fdSLois Curfman McInnes the users manual. 21783e2fdc7SBarry Smith 21836851e7fSLois Curfman McInnes Level: beginner 21936851e7fSLois Curfman McInnes 2209b94acceSBarry Smith .keywords: SNES, nonlinear, set, options, database 2219b94acceSBarry Smith 22215091d37SBarry Smith .seealso: SNESPrintHelp(), SNESSetOptionsPrefix(), SNESSetTypeFromOptions() 2239b94acceSBarry Smith @*/ 2249b94acceSBarry Smith int SNESSetFromOptions(SNES snes) 2259b94acceSBarry Smith { 2269b94acceSBarry Smith double tmp; 2279b94acceSBarry Smith SLES sles; 228f1af5d2fSBarry Smith PetscTruth flg; 229f1af5d2fSBarry Smith int ierr,i,loc[4],nmax = 4; 2303c7409f5SSatish Balay int version = PETSC_DEFAULT; 2319b94acceSBarry Smith double rtol_0 = PETSC_DEFAULT; 2329b94acceSBarry Smith double rtol_max = PETSC_DEFAULT; 2339b94acceSBarry Smith double gamma2 = PETSC_DEFAULT; 2349b94acceSBarry Smith double alpha = PETSC_DEFAULT; 2359b94acceSBarry Smith double alpha2 = PETSC_DEFAULT; 2369b94acceSBarry Smith double threshold = PETSC_DEFAULT; 2379b94acceSBarry Smith 2383a40ed3dSBarry Smith PetscFunctionBegin; 23977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 24015091d37SBarry Smith ierr = SNESSetTypeFromOptions(snes);CHKERRQ(ierr); 241ca161407SBarry Smith 24215091d37SBarry Smith loc[0] = PETSC_DECIDE; loc[1] = PETSC_DECIDE; loc[2] = 300; loc[3] = 300; 2433c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_stol",&tmp, &flg);CHKERRQ(ierr); 244d64ed03dSBarry Smith if (flg) { 245d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 246d64ed03dSBarry Smith } 2473c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_atol",&tmp, &flg);CHKERRQ(ierr); 248d64ed03dSBarry Smith if (flg) { 249d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 250d64ed03dSBarry Smith } 2513c7409f5SSatish Balay ierr = OptionsGetDouble(snes->prefix,"-snes_rtol",&tmp, &flg);CHKERRQ(ierr); 252d64ed03dSBarry Smith if (flg) { 253d64ed03dSBarry Smith ierr = SNESSetTolerances(snes,PETSC_DEFAULT,tmp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 254d64ed03dSBarry Smith } 255f1af5d2fSBarry Smith ierr = OptionsGetInt(snes->prefix,"-snes_max_it",&snes->max_its, PETSC_NULL);CHKERRQ(ierr); 256f1af5d2fSBarry Smith ierr = OptionsGetInt(snes->prefix,"-snes_max_funcs",&snes->max_funcs, PETSC_NULL);CHKERRQ(ierr); 257d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_trtol",&tmp, &flg);CHKERRQ(ierr); 258d64ed03dSBarry Smith if (flg) { ierr = SNESSetTrustRegionTolerance(snes,tmp);CHKERRQ(ierr); } 259d7a720efSLois Curfman McInnes ierr = OptionsGetDouble(snes->prefix,"-snes_fmin",&tmp, &flg);CHKERRQ(ierr); 260d64ed03dSBarry Smith if (flg) { ierr = SNESSetMinimizationFunctionTolerance(snes,tmp);CHKERRQ(ierr);} 2613c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_ksp_ew_conv", &flg);CHKERRQ(ierr); 2623c7409f5SSatish Balay if (flg) { snes->ksp_ewconv = 1; } 263f1af5d2fSBarry Smith ierr = OptionsGetInt(snes->prefix,"-snes_ksp_ew_version",&version,PETSC_NULL);CHKERRQ(ierr); 264f1af5d2fSBarry Smith ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtol0",&rtol_0,PETSC_NULL);CHKERRQ(ierr); 265f1af5d2fSBarry Smith ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_rtolmax",&rtol_max,PETSC_NULL);CHKERRQ(ierr); 266f1af5d2fSBarry Smith ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_gamma",&gamma2,PETSC_NULL);CHKERRQ(ierr); 267f1af5d2fSBarry Smith ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha",&alpha,PETSC_NULL);CHKERRQ(ierr); 268f1af5d2fSBarry Smith ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_alpha2",&alpha2,PETSC_NULL);CHKERRQ(ierr); 269f1af5d2fSBarry Smith ierr = OptionsGetDouble(snes->prefix,"-snes_ksp_ew_threshold",&threshold,PETSC_NULL);CHKERRQ(ierr); 27093c39befSBarry Smith 27193c39befSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_no_convergence_test",&flg);CHKERRQ(ierr); 27293c39befSBarry Smith if (flg) {snes->converged = 0;} 27393c39befSBarry Smith 2749b94acceSBarry Smith ierr = SNES_KSP_SetParametersEW(snes,version,rtol_0,rtol_max,gamma2,alpha, 2759b94acceSBarry Smith alpha2,threshold);CHKERRQ(ierr); 2765bbad29bSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_cancelmonitors",&flg);CHKERRQ(ierr); 2775cd90555SBarry Smith if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 2783c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_monitor",&flg);CHKERRQ(ierr); 279b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);} 2803c7409f5SSatish Balay ierr = OptionsHasName(snes->prefix,"-snes_smonitor",&flg);CHKERRQ(ierr); 281b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);} 2823f1db9ecSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor",&flg);CHKERRQ(ierr); 283b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 284d132466eSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_vecmonitor_update",&flg);CHKERRQ(ierr); 285b8a78c4aSBarry Smith if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitorUpdate,0,0);CHKERRQ(ierr);} 28681f57ec7SBarry Smith ierr = OptionsGetIntArray(snes->prefix,"-snes_xmonitor",loc,&nmax,&flg);CHKERRQ(ierr); 2873c7409f5SSatish Balay if (flg) { 288*c38d4ed2SBarry Smith int rank; 289*c38d4ed2SBarry Smith ierr = MPI_Comm_rank(snes->comm,&rank);CHKERRQ(ierr); 29017699dbbSLois Curfman McInnes if (!rank) { 291*c38d4ed2SBarry Smith ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); 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" 594454a90a3SBarry Smith static int SNESPublish_Petsc(PetscObject obj) 595e24b481bSBarry Smith { 596aa482453SBarry Smith #if defined(PETSC_HAVE_AMS) 597454a90a3SBarry 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 607454a90a3SBarry 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); 612454a90a3SBarry 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 @*/ 864454a90a3SBarry 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 @*/ 1196454a90a3SBarry 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 @*/ 1290454a90a3SBarry 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 { 1372f1af5d2fSBarry Smith int ierr; 1373f1af5d2fSBarry Smith PetscTruth flg; 13743a40ed3dSBarry Smith 13753a40ed3dSBarry Smith PetscFunctionBegin; 137677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 137777c4ece6SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 1378184914b5SBarry Smith PetscCheckSameComm(snes,x); 13798ddd3da0SLois Curfman McInnes snes->vec_sol = snes->vec_sol_always = x; 13809b94acceSBarry Smith 1381c1f60f51SBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_mf_operator", &flg);CHKERRQ(ierr); 1382c1f60f51SBarry Smith /* 1383c1f60f51SBarry Smith This version replaces the user provided Jacobian matrix with a 1384dfa02198SLois Curfman McInnes matrix-free version but still employs the user-provided preconditioner matrix 1385c1f60f51SBarry Smith */ 1386c1f60f51SBarry Smith if (flg) { 1387c1f60f51SBarry Smith Mat J; 13885a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1389c1f60f51SBarry Smith PLogObjectParent(snes,J); 1390c1f60f51SBarry Smith snes->mfshell = J; 1391c1f60f51SBarry Smith snes->jacobian = J; 1392c2aafc4cSSatish Balay if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 1393c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Jacobian routines\n"); 1394d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 1395c1f60f51SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free operator Hessian routines\n"); 1396d4bb536fSBarry Smith } else { 1397a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free operator option"); 1398d4bb536fSBarry Smith } 13995a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1400c1f60f51SBarry Smith } 1401272ac6f2SLois Curfman McInnes ierr = OptionsHasName(snes->prefix,"-snes_mf", &flg);CHKERRQ(ierr); 1402c1f60f51SBarry Smith /* 1403dfa02198SLois Curfman McInnes This version replaces both the user-provided Jacobian and the user- 1404c1f60f51SBarry Smith provided preconditioner matrix with the default matrix free version. 1405c1f60f51SBarry Smith */ 140631615d04SLois Curfman McInnes if (flg) { 1407272ac6f2SLois Curfman McInnes Mat J; 14085a655dc6SBarry Smith ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1409272ac6f2SLois Curfman McInnes PLogObjectParent(snes,J); 1410272ac6f2SLois Curfman McInnes snes->mfshell = J; 141183e56afdSLois Curfman McInnes if (snes->method_class == SNES_NONLINEAR_EQUATIONS) { 141283e56afdSLois Curfman McInnes ierr = SNESSetJacobian(snes,J,J,0,snes->funP);CHKERRQ(ierr); 141394a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Jacobian routines\n"); 1414d64ed03dSBarry Smith } else if (snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION) { 141583e56afdSLois Curfman McInnes ierr = SNESSetHessian(snes,J,J,0,snes->funP);CHKERRQ(ierr); 141694a424c1SBarry Smith PLogInfo(snes,"SNESSetUp: Setting default matrix-free Hessian routines\n"); 1417d4bb536fSBarry Smith } else { 1418a8c6a408SBarry Smith SETERRQ(PETSC_ERR_SUP,0,"Method class doesn't support matrix-free option"); 1419d4bb536fSBarry Smith } 14205a655dc6SBarry Smith ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1421272ac6f2SLois Curfman McInnes } 14229b94acceSBarry Smith if ((snes->method_class == SNES_NONLINEAR_EQUATIONS)) { 14236831982aSBarry Smith PetscTruth iseqtr; 14246831982aSBarry Smith 1425a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 1426a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetFunction() first"); 14275a655dc6SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1428a8c6a408SBarry Smith if (snes->vec_func == snes->vec_sol) { 1429a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_IDN,0,"Solution vector cannot be function vector"); 1430a8c6a408SBarry Smith } 1431a703fe33SLois Curfman McInnes 1432a703fe33SLois Curfman McInnes /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 14336831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,SNESEQTR,&iseqtr);CHKERRQ(ierr); 14346831982aSBarry Smith if (snes->ksp_ewconv && !iseqtr) { 1435a703fe33SLois Curfman McInnes SLES sles; KSP ksp; 1436a703fe33SLois Curfman McInnes ierr = SNESGetSLES(snes,&sles);CHKERRQ(ierr); 1437a703fe33SLois Curfman McInnes ierr = SLESGetKSP(sles,&ksp);CHKERRQ(ierr); 14385a655dc6SBarry Smith ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,(void *)snes);CHKERRQ(ierr); 1439a703fe33SLois Curfman McInnes } 14409b94acceSBarry Smith } else if ((snes->method_class == SNES_UNCONSTRAINED_MINIMIZATION)) { 1441a8c6a408SBarry Smith if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1442a8c6a408SBarry Smith if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetGradient() first"); 1443a8c6a408SBarry Smith if (!snes->computeumfunction) { 1444a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetMinimizationFunction() first"); 1445a8c6a408SBarry Smith } 14465a655dc6SBarry Smith if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,0,"Must call SNESSetHessian()"); 1447d4bb536fSBarry Smith } else { 1448a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Unknown method class"); 1449d4bb536fSBarry Smith } 1450a703fe33SLois Curfman McInnes if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 145182bf6240SBarry Smith snes->setupcalled = 1; 14523a40ed3dSBarry Smith PetscFunctionReturn(0); 14539b94acceSBarry Smith } 14549b94acceSBarry Smith 14555615d1e5SSatish Balay #undef __FUNC__ 1456d4bb536fSBarry Smith #define __FUNC__ "SNESDestroy" 14579b94acceSBarry Smith /*@C 14589b94acceSBarry Smith SNESDestroy - Destroys the nonlinear solver context that was created 14599b94acceSBarry Smith with SNESCreate(). 14609b94acceSBarry Smith 1461c7afd0dbSLois Curfman McInnes Collective on SNES 1462c7afd0dbSLois Curfman McInnes 14639b94acceSBarry Smith Input Parameter: 14649b94acceSBarry Smith . snes - the SNES context 14659b94acceSBarry Smith 146636851e7fSLois Curfman McInnes Level: beginner 146736851e7fSLois Curfman McInnes 14689b94acceSBarry Smith .keywords: SNES, nonlinear, destroy 14699b94acceSBarry Smith 147063a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESSolve() 14719b94acceSBarry Smith @*/ 14729b94acceSBarry Smith int SNESDestroy(SNES snes) 14739b94acceSBarry Smith { 1474b8a78c4aSBarry Smith int i,ierr; 14753a40ed3dSBarry Smith 14763a40ed3dSBarry Smith PetscFunctionBegin; 147777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 14783a40ed3dSBarry Smith if (--snes->refct > 0) PetscFunctionReturn(0); 1479d4bb536fSBarry Smith 1480be0abb6dSBarry Smith /* if memory was published with AMS then destroy it */ 14810f5bd95cSBarry Smith ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1482be0abb6dSBarry Smith 1483e1311b90SBarry Smith if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1484606d414cSSatish Balay if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 1485522c5e43SBarry Smith if (snes->mfshell) {ierr = MatDestroy(snes->mfshell);CHKERRQ(ierr);} 14869b94acceSBarry Smith ierr = SLESDestroy(snes->sles);CHKERRQ(ierr); 1487522c5e43SBarry Smith if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1488b8a78c4aSBarry Smith for (i=0; i<snes->numbermonitors; i++ ) { 1489b8a78c4aSBarry Smith if (snes->monitordestroy[i]) { 1490b8a78c4aSBarry Smith ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1491b8a78c4aSBarry Smith } 1492b8a78c4aSBarry Smith } 14939b94acceSBarry Smith PLogObjectDestroy((PetscObject)snes); 14940452661fSBarry Smith PetscHeaderDestroy((PetscObject)snes); 14953a40ed3dSBarry Smith PetscFunctionReturn(0); 14969b94acceSBarry Smith } 14979b94acceSBarry Smith 14989b94acceSBarry Smith /* ----------- Routines to set solver parameters ---------- */ 14999b94acceSBarry Smith 15005615d1e5SSatish Balay #undef __FUNC__ 15015615d1e5SSatish Balay #define __FUNC__ "SNESSetTolerances" 15029b94acceSBarry Smith /*@ 1503d7a720efSLois Curfman McInnes SNESSetTolerances - Sets various parameters used in convergence tests. 15049b94acceSBarry Smith 1505c7afd0dbSLois Curfman McInnes Collective on SNES 1506c7afd0dbSLois Curfman McInnes 15079b94acceSBarry Smith Input Parameters: 1508c7afd0dbSLois Curfman McInnes + snes - the SNES context 150933174efeSLois Curfman McInnes . atol - absolute convergence tolerance 151033174efeSLois Curfman McInnes . rtol - relative convergence tolerance 151133174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 151233174efeSLois Curfman McInnes of the change in the solution between steps 151333174efeSLois Curfman McInnes . maxit - maximum number of iterations 1514c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1515fee21e36SBarry Smith 151633174efeSLois Curfman McInnes Options Database Keys: 1517c7afd0dbSLois Curfman McInnes + -snes_atol <atol> - Sets atol 1518c7afd0dbSLois Curfman McInnes . -snes_rtol <rtol> - Sets rtol 1519c7afd0dbSLois Curfman McInnes . -snes_stol <stol> - Sets stol 1520c7afd0dbSLois Curfman McInnes . -snes_max_it <maxit> - Sets maxit 1521c7afd0dbSLois Curfman McInnes - -snes_max_funcs <maxf> - Sets maxf 15229b94acceSBarry Smith 1523d7a720efSLois Curfman McInnes Notes: 15249b94acceSBarry Smith The default maximum number of iterations is 50. 15259b94acceSBarry Smith The default maximum number of function evaluations is 1000. 15269b94acceSBarry Smith 152736851e7fSLois Curfman McInnes Level: intermediate 152836851e7fSLois Curfman McInnes 152933174efeSLois Curfman McInnes .keywords: SNES, nonlinear, set, convergence, tolerances 15309b94acceSBarry Smith 1531d7a720efSLois Curfman McInnes .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 15329b94acceSBarry Smith @*/ 1533d7a720efSLois Curfman McInnes int SNESSetTolerances(SNES snes,double atol,double rtol,double stol,int maxit,int maxf) 15349b94acceSBarry Smith { 15353a40ed3dSBarry Smith PetscFunctionBegin; 153677c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1537d7a720efSLois Curfman McInnes if (atol != PETSC_DEFAULT) snes->atol = atol; 1538d7a720efSLois Curfman McInnes if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1539d7a720efSLois Curfman McInnes if (stol != PETSC_DEFAULT) snes->xtol = stol; 1540d7a720efSLois Curfman McInnes if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1541d7a720efSLois Curfman McInnes if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 15423a40ed3dSBarry Smith PetscFunctionReturn(0); 15439b94acceSBarry Smith } 15449b94acceSBarry Smith 15455615d1e5SSatish Balay #undef __FUNC__ 15465615d1e5SSatish Balay #define __FUNC__ "SNESGetTolerances" 15479b94acceSBarry Smith /*@ 154833174efeSLois Curfman McInnes SNESGetTolerances - Gets various parameters used in convergence tests. 154933174efeSLois Curfman McInnes 1550c7afd0dbSLois Curfman McInnes Not Collective 1551c7afd0dbSLois Curfman McInnes 155233174efeSLois Curfman McInnes Input Parameters: 1553c7afd0dbSLois Curfman McInnes + snes - the SNES context 155433174efeSLois Curfman McInnes . atol - absolute convergence tolerance 155533174efeSLois Curfman McInnes . rtol - relative convergence tolerance 155633174efeSLois Curfman McInnes . stol - convergence tolerance in terms of the norm 155733174efeSLois Curfman McInnes of the change in the solution between steps 155833174efeSLois Curfman McInnes . maxit - maximum number of iterations 1559c7afd0dbSLois Curfman McInnes - maxf - maximum number of function evaluations 1560fee21e36SBarry Smith 156133174efeSLois Curfman McInnes Notes: 156233174efeSLois Curfman McInnes The user can specify PETSC_NULL for any parameter that is not needed. 156333174efeSLois Curfman McInnes 156436851e7fSLois Curfman McInnes Level: intermediate 156536851e7fSLois Curfman McInnes 156633174efeSLois Curfman McInnes .keywords: SNES, nonlinear, get, convergence, tolerances 156733174efeSLois Curfman McInnes 156833174efeSLois Curfman McInnes .seealso: SNESSetTolerances() 156933174efeSLois Curfman McInnes @*/ 157033174efeSLois Curfman McInnes int SNESGetTolerances(SNES snes,double *atol,double *rtol,double *stol,int *maxit,int *maxf) 157133174efeSLois Curfman McInnes { 15723a40ed3dSBarry Smith PetscFunctionBegin; 157333174efeSLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 157433174efeSLois Curfman McInnes if (atol) *atol = snes->atol; 157533174efeSLois Curfman McInnes if (rtol) *rtol = snes->rtol; 157633174efeSLois Curfman McInnes if (stol) *stol = snes->xtol; 157733174efeSLois Curfman McInnes if (maxit) *maxit = snes->max_its; 157833174efeSLois Curfman McInnes if (maxf) *maxf = snes->max_funcs; 15793a40ed3dSBarry Smith PetscFunctionReturn(0); 158033174efeSLois Curfman McInnes } 158133174efeSLois Curfman McInnes 15825615d1e5SSatish Balay #undef __FUNC__ 15835615d1e5SSatish Balay #define __FUNC__ "SNESSetTrustRegionTolerance" 158433174efeSLois Curfman McInnes /*@ 15859b94acceSBarry Smith SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 15869b94acceSBarry Smith 1587fee21e36SBarry Smith Collective on SNES 1588fee21e36SBarry Smith 1589c7afd0dbSLois Curfman McInnes Input Parameters: 1590c7afd0dbSLois Curfman McInnes + snes - the SNES context 1591c7afd0dbSLois Curfman McInnes - tol - tolerance 1592c7afd0dbSLois Curfman McInnes 15939b94acceSBarry Smith Options Database Key: 1594c7afd0dbSLois Curfman McInnes . -snes_trtol <tol> - Sets tol 15959b94acceSBarry Smith 159636851e7fSLois Curfman McInnes Level: intermediate 159736851e7fSLois Curfman McInnes 15989b94acceSBarry Smith .keywords: SNES, nonlinear, set, trust region, tolerance 15999b94acceSBarry Smith 1600d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 16019b94acceSBarry Smith @*/ 16029b94acceSBarry Smith int SNESSetTrustRegionTolerance(SNES snes,double tol) 16039b94acceSBarry Smith { 16043a40ed3dSBarry Smith PetscFunctionBegin; 160577c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16069b94acceSBarry Smith snes->deltatol = tol; 16073a40ed3dSBarry Smith PetscFunctionReturn(0); 16089b94acceSBarry Smith } 16099b94acceSBarry Smith 16105615d1e5SSatish Balay #undef __FUNC__ 16115615d1e5SSatish Balay #define __FUNC__ "SNESSetMinimizationFunctionTolerance" 16129b94acceSBarry Smith /*@ 161377c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance - Sets the minimum allowable function tolerance 16149b94acceSBarry Smith for unconstrained minimization solvers. 16159b94acceSBarry Smith 1616fee21e36SBarry Smith Collective on SNES 1617fee21e36SBarry Smith 1618c7afd0dbSLois Curfman McInnes Input Parameters: 1619c7afd0dbSLois Curfman McInnes + snes - the SNES context 1620c7afd0dbSLois Curfman McInnes - ftol - minimum function tolerance 1621c7afd0dbSLois Curfman McInnes 16229b94acceSBarry Smith Options Database Key: 1623c7afd0dbSLois Curfman McInnes . -snes_fmin <ftol> - Sets ftol 16249b94acceSBarry Smith 16259b94acceSBarry Smith Note: 162677c4ece6SBarry Smith SNESSetMinimizationFunctionTolerance() is valid for SNES_UNCONSTRAINED_MINIMIZATION 16279b94acceSBarry Smith methods only. 16289b94acceSBarry Smith 162936851e7fSLois Curfman McInnes Level: intermediate 163036851e7fSLois Curfman McInnes 16319b94acceSBarry Smith .keywords: SNES, nonlinear, set, minimum, convergence, function, tolerance 16329b94acceSBarry Smith 1633d7a720efSLois Curfman McInnes .seealso: SNESSetTolerances(), SNESSetTrustRegionTolerance() 16349b94acceSBarry Smith @*/ 163577c4ece6SBarry Smith int SNESSetMinimizationFunctionTolerance(SNES snes,double ftol) 16369b94acceSBarry Smith { 16373a40ed3dSBarry Smith PetscFunctionBegin; 163877c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 16399b94acceSBarry Smith snes->fmin = ftol; 16403a40ed3dSBarry Smith PetscFunctionReturn(0); 16419b94acceSBarry Smith } 1642df9fa365SBarry Smith /* 1643df9fa365SBarry Smith Duplicate the lg monitors for SNES from KSP; for some reason with 1644df9fa365SBarry Smith dynamic libraries things don't work under Sun4 if we just use 1645df9fa365SBarry Smith macros instead of functions 1646df9fa365SBarry Smith */ 1647ce1608b8SBarry Smith #undef __FUNC__ 1648ce1608b8SBarry Smith #define __FUNC__ "SNESLGMonitor" 1649ce1608b8SBarry Smith int SNESLGMonitor(SNES snes,int it,double norm,void *ctx) 1650ce1608b8SBarry Smith { 1651ce1608b8SBarry Smith int ierr; 1652ce1608b8SBarry Smith 1653ce1608b8SBarry Smith PetscFunctionBegin; 1654184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1655ce1608b8SBarry Smith ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1656ce1608b8SBarry Smith PetscFunctionReturn(0); 1657ce1608b8SBarry Smith } 1658ce1608b8SBarry Smith 1659df9fa365SBarry Smith #undef __FUNC__ 1660df9fa365SBarry Smith #define __FUNC__ "SNESLGMonitorCreate" 1661df9fa365SBarry Smith int SNESLGMonitorCreate(char *host,char *label,int x,int y,int m,int n, DrawLG *draw) 1662df9fa365SBarry Smith { 1663df9fa365SBarry Smith int ierr; 1664df9fa365SBarry Smith 1665df9fa365SBarry Smith PetscFunctionBegin; 1666df9fa365SBarry Smith ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1667df9fa365SBarry Smith PetscFunctionReturn(0); 1668df9fa365SBarry Smith } 1669df9fa365SBarry Smith 1670df9fa365SBarry Smith #undef __FUNC__ 1671df9fa365SBarry Smith #define __FUNC__ "SNESLGMonitorDestroy" 1672df9fa365SBarry Smith int SNESLGMonitorDestroy(DrawLG draw) 1673df9fa365SBarry Smith { 1674df9fa365SBarry Smith int ierr; 1675df9fa365SBarry Smith 1676df9fa365SBarry Smith PetscFunctionBegin; 1677df9fa365SBarry Smith ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1678df9fa365SBarry Smith PetscFunctionReturn(0); 1679df9fa365SBarry Smith } 1680df9fa365SBarry Smith 16819b94acceSBarry Smith /* ------------ Routines to set performance monitoring options ----------- */ 16829b94acceSBarry Smith 16835615d1e5SSatish Balay #undef __FUNC__ 1684d4bb536fSBarry Smith #define __FUNC__ "SNESSetMonitor" 16859b94acceSBarry Smith /*@C 16865cd90555SBarry Smith SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 16879b94acceSBarry Smith iteration of the nonlinear solver to display the iteration's 16889b94acceSBarry Smith progress. 16899b94acceSBarry Smith 1690fee21e36SBarry Smith Collective on SNES 1691fee21e36SBarry Smith 1692c7afd0dbSLois Curfman McInnes Input Parameters: 1693c7afd0dbSLois Curfman McInnes + snes - the SNES context 1694c7afd0dbSLois Curfman McInnes . func - monitoring routine 1695b8a78c4aSBarry Smith . mctx - [optional] user-defined context for private data for the 1696c7afd0dbSLois Curfman McInnes monitor routine (may be PETSC_NULL) 1697b8a78c4aSBarry Smith - monitordestroy - options routine that frees monitor context 16989b94acceSBarry Smith 1699c7afd0dbSLois Curfman McInnes Calling sequence of func: 170040a0c1c6SLois Curfman McInnes $ int func(SNES snes,int its, double norm,void *mctx) 1701c7afd0dbSLois Curfman McInnes 1702c7afd0dbSLois Curfman McInnes + snes - the SNES context 1703c7afd0dbSLois Curfman McInnes . its - iteration number 1704c7afd0dbSLois Curfman McInnes . norm - 2-norm function value (may be estimated) 1705c7afd0dbSLois Curfman McInnes for SNES_NONLINEAR_EQUATIONS methods 170640a0c1c6SLois Curfman McInnes . norm - 2-norm gradient value (may be estimated) 1707c7afd0dbSLois Curfman McInnes for SNES_UNCONSTRAINED_MINIMIZATION methods 170840a0c1c6SLois Curfman McInnes - mctx - [optional] monitoring context 17099b94acceSBarry Smith 17109665c990SLois Curfman McInnes Options Database Keys: 1711c7afd0dbSLois Curfman McInnes + -snes_monitor - sets SNESDefaultMonitor() 1712c7afd0dbSLois Curfman McInnes . -snes_xmonitor - sets line graph monitor, 1713c7afd0dbSLois Curfman McInnes uses SNESLGMonitorCreate() 1714c7afd0dbSLois Curfman McInnes _ -snes_cancelmonitors - cancels all monitors that have 1715c7afd0dbSLois Curfman McInnes been hardwired into a code by 1716c7afd0dbSLois Curfman McInnes calls to SNESSetMonitor(), but 1717c7afd0dbSLois Curfman McInnes does not cancel those set via 1718c7afd0dbSLois Curfman McInnes the options database. 17199665c990SLois Curfman McInnes 1720639f9d9dSBarry Smith Notes: 17216bc08f3fSLois Curfman McInnes Several different monitoring routines may be set by calling 17226bc08f3fSLois Curfman McInnes SNESSetMonitor() multiple times; all will be called in the 17236bc08f3fSLois Curfman McInnes order in which they were set. 1724639f9d9dSBarry Smith 172536851e7fSLois Curfman McInnes Level: intermediate 172636851e7fSLois Curfman McInnes 17279b94acceSBarry Smith .keywords: SNES, nonlinear, set, monitor 17289b94acceSBarry Smith 17295cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESClearMonitor() 17309b94acceSBarry Smith @*/ 1731b8a78c4aSBarry Smith int SNESSetMonitor( SNES snes, int (*func)(SNES,int,double,void*),void *mctx,int (*monitordestroy)(void *)) 17329b94acceSBarry Smith { 17333a40ed3dSBarry Smith PetscFunctionBegin; 1734184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1735639f9d9dSBarry Smith if (snes->numbermonitors >= MAXSNESMONITORS) { 1736a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,0,"Too many monitors set"); 1737639f9d9dSBarry Smith } 1738639f9d9dSBarry Smith 1739639f9d9dSBarry Smith snes->monitor[snes->numbermonitors] = func; 1740b8a78c4aSBarry Smith snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1741639f9d9dSBarry Smith snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 17423a40ed3dSBarry Smith PetscFunctionReturn(0); 17439b94acceSBarry Smith } 17449b94acceSBarry Smith 17455615d1e5SSatish Balay #undef __FUNC__ 17465cd90555SBarry Smith #define __FUNC__ "SNESClearMonitor" 17475cd90555SBarry Smith /*@C 17485cd90555SBarry Smith SNESClearMonitor - Clears all the monitor functions for a SNES object. 17495cd90555SBarry Smith 1750c7afd0dbSLois Curfman McInnes Collective on SNES 1751c7afd0dbSLois Curfman McInnes 17525cd90555SBarry Smith Input Parameters: 17535cd90555SBarry Smith . snes - the SNES context 17545cd90555SBarry Smith 17555cd90555SBarry Smith Options Database: 1756c7afd0dbSLois Curfman McInnes . -snes_cancelmonitors - cancels all monitors that have been hardwired 1757c7afd0dbSLois Curfman McInnes into a code by calls to SNESSetMonitor(), but does not cancel those 1758c7afd0dbSLois Curfman McInnes set via the options database 17595cd90555SBarry Smith 17605cd90555SBarry Smith Notes: 17615cd90555SBarry Smith There is no way to clear one specific monitor from a SNES object. 17625cd90555SBarry Smith 176336851e7fSLois Curfman McInnes Level: intermediate 176436851e7fSLois Curfman McInnes 17655cd90555SBarry Smith .keywords: SNES, nonlinear, set, monitor 17665cd90555SBarry Smith 17675cd90555SBarry Smith .seealso: SNESDefaultMonitor(), SNESSetMonitor() 17685cd90555SBarry Smith @*/ 17695cd90555SBarry Smith int SNESClearMonitor( SNES snes ) 17705cd90555SBarry Smith { 17715cd90555SBarry Smith PetscFunctionBegin; 1772184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 17735cd90555SBarry Smith snes->numbermonitors = 0; 17745cd90555SBarry Smith PetscFunctionReturn(0); 17755cd90555SBarry Smith } 17765cd90555SBarry Smith 17775cd90555SBarry Smith #undef __FUNC__ 1778d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceTest" 17799b94acceSBarry Smith /*@C 17809b94acceSBarry Smith SNESSetConvergenceTest - Sets the function that is to be used 17819b94acceSBarry Smith to test for convergence of the nonlinear iterative solution. 17829b94acceSBarry Smith 1783fee21e36SBarry Smith Collective on SNES 1784fee21e36SBarry Smith 1785c7afd0dbSLois Curfman McInnes Input Parameters: 1786c7afd0dbSLois Curfman McInnes + snes - the SNES context 1787c7afd0dbSLois Curfman McInnes . func - routine to test for convergence 1788c7afd0dbSLois Curfman McInnes - cctx - [optional] context for private data for the convergence routine 1789c7afd0dbSLois Curfman McInnes (may be PETSC_NULL) 17909b94acceSBarry Smith 1791c7afd0dbSLois Curfman McInnes Calling sequence of func: 1792184914b5SBarry Smith $ int func (SNES snes,double xnorm,double gnorm,double f,SNESConvergedReason *reason,void *cctx) 1793c7afd0dbSLois Curfman McInnes 1794c7afd0dbSLois Curfman McInnes + snes - the SNES context 1795c7afd0dbSLois Curfman McInnes . cctx - [optional] convergence context 1796184914b5SBarry Smith . reason - reason for convergence/divergence 1797c7afd0dbSLois Curfman McInnes . xnorm - 2-norm of current iterate 1798c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current step (SNES_NONLINEAR_EQUATIONS methods) 1799c7afd0dbSLois Curfman McInnes . f - 2-norm of function (SNES_NONLINEAR_EQUATIONS methods) 1800c7afd0dbSLois Curfman McInnes . gnorm - 2-norm of current gradient (SNES_UNCONSTRAINED_MINIMIZATION methods) 1801c7afd0dbSLois Curfman McInnes - f - function value (SNES_UNCONSTRAINED_MINIMIZATION methods) 18029b94acceSBarry Smith 180336851e7fSLois Curfman McInnes Level: advanced 180436851e7fSLois Curfman McInnes 18059b94acceSBarry Smith .keywords: SNES, nonlinear, set, convergence, test 18069b94acceSBarry Smith 180740191667SLois Curfman McInnes .seealso: SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 180840191667SLois Curfman McInnes SNESConverged_UM_LS(), SNESConverged_UM_TR() 18099b94acceSBarry Smith @*/ 1810184914b5SBarry Smith int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,double,double,double,SNESConvergedReason*,void*),void *cctx) 18119b94acceSBarry Smith { 18123a40ed3dSBarry Smith PetscFunctionBegin; 1813184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 18149b94acceSBarry Smith (snes)->converged = func; 18159b94acceSBarry Smith (snes)->cnvP = cctx; 18163a40ed3dSBarry Smith PetscFunctionReturn(0); 18179b94acceSBarry Smith } 18189b94acceSBarry Smith 18195615d1e5SSatish Balay #undef __FUNC__ 1820184914b5SBarry Smith #define __FUNC__ "SNESGetConvergedReason" 1821184914b5SBarry Smith /*@C 1822184914b5SBarry Smith SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1823184914b5SBarry Smith 1824184914b5SBarry Smith Not Collective 1825184914b5SBarry Smith 1826184914b5SBarry Smith Input Parameter: 1827184914b5SBarry Smith . snes - the SNES context 1828184914b5SBarry Smith 1829184914b5SBarry Smith Output Parameter: 1830184914b5SBarry Smith . reason - negative value indicates diverged, positive value converged, see snes.h or the 1831184914b5SBarry Smith manual pages for the individual convergence tests for complete lists 1832184914b5SBarry Smith 1833184914b5SBarry Smith Level: intermediate 1834184914b5SBarry Smith 1835184914b5SBarry Smith Notes: Can only be called after the call the SNESSolve() is complete. 1836184914b5SBarry Smith 1837184914b5SBarry Smith .keywords: SNES, nonlinear, set, convergence, test 1838184914b5SBarry Smith 1839184914b5SBarry Smith .seealso: SNESSetConvergenceTest(), SNESConverged_EQ_LS(), SNESConverged_EQ_TR(), 1840184914b5SBarry Smith SNESConverged_UM_LS(), SNESConverged_UM_TR() 1841184914b5SBarry Smith @*/ 1842184914b5SBarry Smith int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1843184914b5SBarry Smith { 1844184914b5SBarry Smith PetscFunctionBegin; 1845184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1846184914b5SBarry Smith *reason = snes->reason; 1847184914b5SBarry Smith PetscFunctionReturn(0); 1848184914b5SBarry Smith } 1849184914b5SBarry Smith 1850184914b5SBarry Smith #undef __FUNC__ 1851d4bb536fSBarry Smith #define __FUNC__ "SNESSetConvergenceHistory" 1852c9005455SLois Curfman McInnes /*@ 1853c9005455SLois Curfman McInnes SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1854c9005455SLois Curfman McInnes 1855fee21e36SBarry Smith Collective on SNES 1856fee21e36SBarry Smith 1857c7afd0dbSLois Curfman McInnes Input Parameters: 1858c7afd0dbSLois Curfman McInnes + snes - iterative context obtained from SNESCreate() 1859c7afd0dbSLois Curfman McInnes . a - array to hold history 1860758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1861758f92a0SBarry Smith negative if not converged) for each solve. 1862758f92a0SBarry Smith . na - size of a and its 1863758f92a0SBarry Smith - reset - PETSC_TRUTH indicates each new nonlinear solve resets the history counter to zero, 1864758f92a0SBarry Smith else it continues storing new values for new nonlinear solves after the old ones 1865c7afd0dbSLois Curfman McInnes 1866c9005455SLois Curfman McInnes Notes: 1867c9005455SLois Curfman McInnes If set, this array will contain the function norms (for 1868c9005455SLois Curfman McInnes SNES_NONLINEAR_EQUATIONS methods) or gradient norms 1869c9005455SLois Curfman McInnes (for SNES_UNCONSTRAINED_MINIMIZATION methods) computed 1870c9005455SLois Curfman McInnes at each step. 1871c9005455SLois Curfman McInnes 1872c9005455SLois Curfman McInnes This routine is useful, e.g., when running a code for purposes 1873c9005455SLois Curfman McInnes of accurate performance monitoring, when no I/O should be done 1874c9005455SLois Curfman McInnes during the section of code that is being timed. 1875c9005455SLois Curfman McInnes 187636851e7fSLois Curfman McInnes Level: intermediate 187736851e7fSLois Curfman McInnes 1878c9005455SLois Curfman McInnes .keywords: SNES, set, convergence, history 1879758f92a0SBarry Smith 188008405cd6SLois Curfman McInnes .seealso: SNESGetConvergenceHistory() 1881758f92a0SBarry Smith 1882c9005455SLois Curfman McInnes @*/ 1883758f92a0SBarry Smith int SNESSetConvergenceHistory(SNES snes, double *a, int *its,int na,PetscTruth reset) 1884c9005455SLois Curfman McInnes { 18853a40ed3dSBarry Smith PetscFunctionBegin; 1886c9005455SLois Curfman McInnes PetscValidHeaderSpecific(snes,SNES_COOKIE); 1887c9005455SLois Curfman McInnes if (na) PetscValidScalarPointer(a); 1888c9005455SLois Curfman McInnes snes->conv_hist = a; 1889758f92a0SBarry Smith snes->conv_hist_its = its; 1890758f92a0SBarry Smith snes->conv_hist_max = na; 1891758f92a0SBarry Smith snes->conv_hist_reset = reset; 1892758f92a0SBarry Smith PetscFunctionReturn(0); 1893758f92a0SBarry Smith } 1894758f92a0SBarry Smith 1895758f92a0SBarry Smith #undef __FUNC__ 1896758f92a0SBarry Smith #define __FUNC__ "SNESGetConvergenceHistory" 18970c4c9dddSBarry Smith /*@C 1898758f92a0SBarry Smith SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1899758f92a0SBarry Smith 1900758f92a0SBarry Smith Collective on SNES 1901758f92a0SBarry Smith 1902758f92a0SBarry Smith Input Parameter: 1903758f92a0SBarry Smith . snes - iterative context obtained from SNESCreate() 1904758f92a0SBarry Smith 1905758f92a0SBarry Smith Output Parameters: 1906758f92a0SBarry Smith . a - array to hold history 1907758f92a0SBarry Smith . its - integer array holds the number of linear iterations (or 1908758f92a0SBarry Smith negative if not converged) for each solve. 1909758f92a0SBarry Smith - na - size of a and its 1910758f92a0SBarry Smith 1911758f92a0SBarry Smith Notes: 1912758f92a0SBarry Smith The calling sequence for this routine in Fortran is 1913758f92a0SBarry Smith $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1914758f92a0SBarry Smith 1915758f92a0SBarry Smith This routine is useful, e.g., when running a code for purposes 1916758f92a0SBarry Smith of accurate performance monitoring, when no I/O should be done 1917758f92a0SBarry Smith during the section of code that is being timed. 1918758f92a0SBarry Smith 1919758f92a0SBarry Smith Level: intermediate 1920758f92a0SBarry Smith 1921758f92a0SBarry Smith .keywords: SNES, get, convergence, history 1922758f92a0SBarry Smith 1923758f92a0SBarry Smith .seealso: SNESSetConvergencHistory() 1924758f92a0SBarry Smith 1925758f92a0SBarry Smith @*/ 1926758f92a0SBarry Smith int SNESGetConvergenceHistory(SNES snes, double **a, int **its,int *na) 1927758f92a0SBarry Smith { 1928758f92a0SBarry Smith PetscFunctionBegin; 1929758f92a0SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1930758f92a0SBarry Smith if (a) *a = snes->conv_hist; 1931758f92a0SBarry Smith if (its) *its = snes->conv_hist_its; 1932758f92a0SBarry Smith if (na) *na = snes->conv_hist_len; 19333a40ed3dSBarry Smith PetscFunctionReturn(0); 1934c9005455SLois Curfman McInnes } 1935c9005455SLois Curfman McInnes 1936c9005455SLois Curfman McInnes #undef __FUNC__ 19375615d1e5SSatish Balay #define __FUNC__ "SNESScaleStep_Private" 19389b94acceSBarry Smith /* 19399b94acceSBarry Smith SNESScaleStep_Private - Scales a step so that its length is less than the 19409b94acceSBarry Smith positive parameter delta. 19419b94acceSBarry Smith 19429b94acceSBarry Smith Input Parameters: 1943c7afd0dbSLois Curfman McInnes + snes - the SNES context 19449b94acceSBarry Smith . y - approximate solution of linear system 19459b94acceSBarry Smith . fnorm - 2-norm of current function 1946c7afd0dbSLois Curfman McInnes - delta - trust region size 19479b94acceSBarry Smith 19489b94acceSBarry Smith Output Parameters: 1949c7afd0dbSLois Curfman McInnes + gpnorm - predicted function norm at the new point, assuming local 19509b94acceSBarry Smith linearization. The value is zero if the step lies within the trust 19519b94acceSBarry Smith region, and exceeds zero otherwise. 1952c7afd0dbSLois Curfman McInnes - ynorm - 2-norm of the step 19539b94acceSBarry Smith 19549b94acceSBarry Smith Note: 19556831982aSBarry Smith For non-trust region methods such as SNESEQLS, the parameter delta 19569b94acceSBarry Smith is set to be the maximum allowable step size. 19579b94acceSBarry Smith 19589b94acceSBarry Smith .keywords: SNES, nonlinear, scale, step 19599b94acceSBarry Smith */ 19609b94acceSBarry Smith int SNESScaleStep_Private(SNES snes,Vec y,double *fnorm,double *delta, 19619b94acceSBarry Smith double *gpnorm,double *ynorm) 19629b94acceSBarry Smith { 19639b94acceSBarry Smith double norm; 19649b94acceSBarry Smith Scalar cnorm; 19653a40ed3dSBarry Smith int ierr; 19663a40ed3dSBarry Smith 19673a40ed3dSBarry Smith PetscFunctionBegin; 1968184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 1969184914b5SBarry Smith PetscValidHeaderSpecific(y,VEC_COOKIE); 1970184914b5SBarry Smith PetscCheckSameComm(snes,y); 1971184914b5SBarry Smith 19723a40ed3dSBarry Smith ierr = VecNorm(y,NORM_2, &norm );CHKERRQ(ierr); 19739b94acceSBarry Smith if (norm > *delta) { 19749b94acceSBarry Smith norm = *delta/norm; 19759b94acceSBarry Smith *gpnorm = (1.0 - norm)*(*fnorm); 19769b94acceSBarry Smith cnorm = norm; 19779b94acceSBarry Smith VecScale( &cnorm, y ); 19789b94acceSBarry Smith *ynorm = *delta; 19799b94acceSBarry Smith } else { 19809b94acceSBarry Smith *gpnorm = 0.0; 19819b94acceSBarry Smith *ynorm = norm; 19829b94acceSBarry Smith } 19833a40ed3dSBarry Smith PetscFunctionReturn(0); 19849b94acceSBarry Smith } 19859b94acceSBarry Smith 19865615d1e5SSatish Balay #undef __FUNC__ 19875615d1e5SSatish Balay #define __FUNC__ "SNESSolve" 19889b94acceSBarry Smith /*@ 19899b94acceSBarry Smith SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 199063a78c88SLois Curfman McInnes SNESCreate() and optional routines of the form SNESSetXXX(). 19919b94acceSBarry Smith 1992c7afd0dbSLois Curfman McInnes Collective on SNES 1993c7afd0dbSLois Curfman McInnes 1994b2002411SLois Curfman McInnes Input Parameters: 1995c7afd0dbSLois Curfman McInnes + snes - the SNES context 1996c7afd0dbSLois Curfman McInnes - x - the solution vector 19979b94acceSBarry Smith 19989b94acceSBarry Smith Output Parameter: 1999b2002411SLois Curfman McInnes . its - number of iterations until termination 20009b94acceSBarry Smith 2001b2002411SLois Curfman McInnes Notes: 20028ddd3da0SLois Curfman McInnes The user should initialize the vector, x, with the initial guess 20038ddd3da0SLois Curfman McInnes for the nonlinear solve prior to calling SNESSolve. In particular, 20048ddd3da0SLois Curfman McInnes to employ an initial guess of zero, the user should explicitly set 20058ddd3da0SLois Curfman McInnes this vector to zero by calling VecSet(). 20068ddd3da0SLois Curfman McInnes 200736851e7fSLois Curfman McInnes Level: beginner 200836851e7fSLois Curfman McInnes 20099b94acceSBarry Smith .keywords: SNES, nonlinear, solve 20109b94acceSBarry Smith 201163a78c88SLois Curfman McInnes .seealso: SNESCreate(), SNESDestroy() 20129b94acceSBarry Smith @*/ 20138ddd3da0SLois Curfman McInnes int SNESSolve(SNES snes,Vec x,int *its) 20149b94acceSBarry Smith { 2015f1af5d2fSBarry Smith int ierr; 2016f1af5d2fSBarry Smith PetscTruth flg; 2017052efed2SBarry Smith 20183a40ed3dSBarry Smith PetscFunctionBegin; 201977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2020184914b5SBarry Smith PetscValidHeaderSpecific(x,VEC_COOKIE); 2021184914b5SBarry Smith PetscCheckSameComm(snes,x); 202274679c65SBarry Smith PetscValidIntPointer(its); 202382bf6240SBarry Smith if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 2024c4fc05e7SBarry Smith else {snes->vec_sol = snes->vec_sol_always = x;} 2025758f92a0SBarry Smith if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 20269b94acceSBarry Smith PLogEventBegin(SNES_Solve,snes,0,0,0); 2027c96a6f78SLois Curfman McInnes snes->nfuncs = 0; snes->linear_its = 0; snes->nfailures = 0; 20289b94acceSBarry Smith ierr = (*(snes)->solve)(snes,its);CHKERRQ(ierr); 20299b94acceSBarry Smith PLogEventEnd(SNES_Solve,snes,0,0,0); 20300f5bd95cSBarry Smith ierr = OptionsHasName(snes->prefix,"-snes_view", &flg);CHKERRQ(ierr); 20316d4a8577SBarry Smith if (flg) { ierr = SNESView(snes,VIEWER_STDOUT_WORLD);CHKERRQ(ierr); } 20323a40ed3dSBarry Smith PetscFunctionReturn(0); 20339b94acceSBarry Smith } 20349b94acceSBarry Smith 20359b94acceSBarry Smith /* --------- Internal routines for SNES Package --------- */ 20369b94acceSBarry Smith 20375615d1e5SSatish Balay #undef __FUNC__ 20385615d1e5SSatish Balay #define __FUNC__ "SNESSetType" 203982bf6240SBarry Smith /*@C 20404b0e389bSBarry Smith SNESSetType - Sets the method for the nonlinear solver. 20419b94acceSBarry Smith 2042fee21e36SBarry Smith Collective on SNES 2043fee21e36SBarry Smith 2044c7afd0dbSLois Curfman McInnes Input Parameters: 2045c7afd0dbSLois Curfman McInnes + snes - the SNES context 2046454a90a3SBarry Smith - type - a known method 2047c7afd0dbSLois Curfman McInnes 2048c7afd0dbSLois Curfman McInnes Options Database Key: 2049454a90a3SBarry Smith . -snes_type <type> - Sets the method; use -help for a list 2050c7afd0dbSLois Curfman McInnes of available methods (for instance, ls or tr) 2051ae12b187SLois Curfman McInnes 20529b94acceSBarry Smith Notes: 20539b94acceSBarry Smith See "petsc/include/snes.h" for available methods (for instance) 20546831982aSBarry Smith + SNESEQLS - Newton's method with line search 2055c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 20566831982aSBarry Smith . SNESEQTR - Newton's method with trust region 2057c7afd0dbSLois Curfman McInnes (systems of nonlinear equations) 20586831982aSBarry Smith . SNESUMTR - Newton's method with trust region 2059c7afd0dbSLois Curfman McInnes (unconstrained minimization) 20606831982aSBarry Smith - SNESUMLS - Newton's method with line search 2061c7afd0dbSLois Curfman McInnes (unconstrained minimization) 20629b94acceSBarry Smith 2063ae12b187SLois Curfman McInnes Normally, it is best to use the SNESSetFromOptions() command and then 2064ae12b187SLois Curfman McInnes set the SNES solver type from the options database rather than by using 2065ae12b187SLois Curfman McInnes this routine. Using the options database provides the user with 2066ae12b187SLois Curfman McInnes maximum flexibility in evaluating the many nonlinear solvers. 2067ae12b187SLois Curfman McInnes The SNESSetType() routine is provided for those situations where it 2068ae12b187SLois Curfman McInnes is necessary to set the nonlinear solver independently of the command 2069ae12b187SLois Curfman McInnes line or options database. This might be the case, for example, when 2070ae12b187SLois Curfman McInnes the choice of solver changes during the execution of the program, 2071ae12b187SLois Curfman McInnes and the user's application is taking responsibility for choosing the 207236851e7fSLois Curfman McInnes appropriate method. In other words, this routine is not for beginners. 207336851e7fSLois Curfman McInnes 207436851e7fSLois Curfman McInnes Level: intermediate 2075a703fe33SLois Curfman McInnes 2076454a90a3SBarry Smith .keywords: SNES, set, type 20779b94acceSBarry Smith @*/ 2078454a90a3SBarry Smith int SNESSetType(SNES snes,SNESType type) 20799b94acceSBarry Smith { 20800f5bd95cSBarry Smith int ierr, (*r)(SNES); 20816831982aSBarry Smith PetscTruth match; 20823a40ed3dSBarry Smith 20833a40ed3dSBarry Smith PetscFunctionBegin; 208477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 20850f5bd95cSBarry Smith PetscValidCharPointer(type); 208682bf6240SBarry Smith 20876831982aSBarry Smith ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 20880f5bd95cSBarry Smith if (match) PetscFunctionReturn(0); 208982bf6240SBarry Smith 209082bf6240SBarry Smith if (snes->setupcalled) { 2091e1311b90SBarry Smith ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 209282bf6240SBarry Smith snes->data = 0; 209382bf6240SBarry Smith } 20947d1a2b2bSBarry Smith 20959b94acceSBarry Smith /* Get the function pointers for the iterative method requested */ 209682bf6240SBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 209782bf6240SBarry Smith 2098454a90a3SBarry Smith ierr = FListFind(snes->comm, SNESList, type,(int (**)(void *)) &r );CHKERRQ(ierr); 209982bf6240SBarry Smith 2100454a90a3SBarry Smith if (!r) SETERRQ1(1,1,"Unable to find requested SNES type %s",type); 21011548b14aSBarry Smith 2102606d414cSSatish Balay if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 210382bf6240SBarry Smith snes->data = 0; 21043a40ed3dSBarry Smith ierr = (*r)(snes);CHKERRQ(ierr); 210582bf6240SBarry Smith 2106454a90a3SBarry Smith ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 210782bf6240SBarry Smith snes->set_method_called = 1; 210882bf6240SBarry Smith 21093a40ed3dSBarry Smith PetscFunctionReturn(0); 21109b94acceSBarry Smith } 21119b94acceSBarry Smith 2112a847f771SSatish Balay 21139b94acceSBarry Smith /* --------------------------------------------------------------------- */ 21145615d1e5SSatish Balay #undef __FUNC__ 2115d4bb536fSBarry Smith #define __FUNC__ "SNESRegisterDestroy" 21169b94acceSBarry Smith /*@C 21179b94acceSBarry Smith SNESRegisterDestroy - Frees the list of nonlinear solvers that were 2118f1af5d2fSBarry Smith registered by SNESRegisterDynamic(). 21199b94acceSBarry Smith 2120fee21e36SBarry Smith Not Collective 2121fee21e36SBarry Smith 212236851e7fSLois Curfman McInnes Level: advanced 212336851e7fSLois Curfman McInnes 21249b94acceSBarry Smith .keywords: SNES, nonlinear, register, destroy 21259b94acceSBarry Smith 21269b94acceSBarry Smith .seealso: SNESRegisterAll(), SNESRegisterAll() 21279b94acceSBarry Smith @*/ 2128cf256101SBarry Smith int SNESRegisterDestroy(void) 21299b94acceSBarry Smith { 213082bf6240SBarry Smith int ierr; 213182bf6240SBarry Smith 21323a40ed3dSBarry Smith PetscFunctionBegin; 213382bf6240SBarry Smith if (SNESList) { 2134488ecbafSBarry Smith ierr = FListDestroy( SNESList );CHKERRQ(ierr); 213582bf6240SBarry Smith SNESList = 0; 21369b94acceSBarry Smith } 213784cb2905SBarry Smith SNESRegisterAllCalled = 0; 21383a40ed3dSBarry Smith PetscFunctionReturn(0); 21399b94acceSBarry Smith } 21409b94acceSBarry Smith 21415615d1e5SSatish Balay #undef __FUNC__ 2142d4bb536fSBarry Smith #define __FUNC__ "SNESGetType" 21439b94acceSBarry Smith /*@C 21449a28b0a6SLois Curfman McInnes SNESGetType - Gets the SNES method type and name (as a string). 21459b94acceSBarry Smith 2146c7afd0dbSLois Curfman McInnes Not Collective 2147c7afd0dbSLois Curfman McInnes 21489b94acceSBarry Smith Input Parameter: 21494b0e389bSBarry Smith . snes - nonlinear solver context 21509b94acceSBarry Smith 21519b94acceSBarry Smith Output Parameter: 2152454a90a3SBarry Smith . type - SNES method (a charactor string) 21539b94acceSBarry Smith 215436851e7fSLois Curfman McInnes Level: intermediate 215536851e7fSLois Curfman McInnes 2156454a90a3SBarry Smith .keywords: SNES, nonlinear, get, type, name 21579b94acceSBarry Smith @*/ 2158454a90a3SBarry Smith int SNESGetType(SNES snes, SNESType *type) 21599b94acceSBarry Smith { 21603a40ed3dSBarry Smith PetscFunctionBegin; 2161184914b5SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2162454a90a3SBarry Smith *type = snes->type_name; 21633a40ed3dSBarry Smith PetscFunctionReturn(0); 21649b94acceSBarry Smith } 21659b94acceSBarry Smith 21665615d1e5SSatish Balay #undef __FUNC__ 2167d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolution" 21689b94acceSBarry Smith /*@C 21699b94acceSBarry Smith SNESGetSolution - Returns the vector where the approximate solution is 21709b94acceSBarry Smith stored. 21719b94acceSBarry Smith 2172c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2173c7afd0dbSLois Curfman McInnes 21749b94acceSBarry Smith Input Parameter: 21759b94acceSBarry Smith . snes - the SNES context 21769b94acceSBarry Smith 21779b94acceSBarry Smith Output Parameter: 21789b94acceSBarry Smith . x - the solution 21799b94acceSBarry Smith 218036851e7fSLois Curfman McInnes Level: advanced 218136851e7fSLois Curfman McInnes 21829b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution 21839b94acceSBarry Smith 2184a86d99e1SLois Curfman McInnes .seealso: SNESGetFunction(), SNESGetGradient(), SNESGetSolutionUpdate() 21859b94acceSBarry Smith @*/ 21869b94acceSBarry Smith int SNESGetSolution(SNES snes,Vec *x) 21879b94acceSBarry Smith { 21883a40ed3dSBarry Smith PetscFunctionBegin; 218977c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 21909b94acceSBarry Smith *x = snes->vec_sol_always; 21913a40ed3dSBarry Smith PetscFunctionReturn(0); 21929b94acceSBarry Smith } 21939b94acceSBarry Smith 21945615d1e5SSatish Balay #undef __FUNC__ 2195d4bb536fSBarry Smith #define __FUNC__ "SNESGetSolutionUpdate" 21969b94acceSBarry Smith /*@C 21979b94acceSBarry Smith SNESGetSolutionUpdate - Returns the vector where the solution update is 21989b94acceSBarry Smith stored. 21999b94acceSBarry Smith 2200c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2201c7afd0dbSLois Curfman McInnes 22029b94acceSBarry Smith Input Parameter: 22039b94acceSBarry Smith . snes - the SNES context 22049b94acceSBarry Smith 22059b94acceSBarry Smith Output Parameter: 22069b94acceSBarry Smith . x - the solution update 22079b94acceSBarry Smith 220836851e7fSLois Curfman McInnes Level: advanced 220936851e7fSLois Curfman McInnes 22109b94acceSBarry Smith .keywords: SNES, nonlinear, get, solution, update 22119b94acceSBarry Smith 22129b94acceSBarry Smith .seealso: SNESGetSolution(), SNESGetFunction 22139b94acceSBarry Smith @*/ 22149b94acceSBarry Smith int SNESGetSolutionUpdate(SNES snes,Vec *x) 22159b94acceSBarry Smith { 22163a40ed3dSBarry Smith PetscFunctionBegin; 221777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 22189b94acceSBarry Smith *x = snes->vec_sol_update_always; 22193a40ed3dSBarry Smith PetscFunctionReturn(0); 22209b94acceSBarry Smith } 22219b94acceSBarry Smith 22225615d1e5SSatish Balay #undef __FUNC__ 2223d4bb536fSBarry Smith #define __FUNC__ "SNESGetFunction" 22249b94acceSBarry Smith /*@C 22253638b69dSLois Curfman McInnes SNESGetFunction - Returns the vector where the function is stored. 22269b94acceSBarry Smith 2227c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2228c7afd0dbSLois Curfman McInnes 22299b94acceSBarry Smith Input Parameter: 22309b94acceSBarry Smith . snes - the SNES context 22319b94acceSBarry Smith 22329b94acceSBarry Smith Output Parameter: 22337bf4e008SBarry Smith + r - the function (or PETSC_NULL) 22347bf4e008SBarry Smith - ctx - the function context (or PETSC_NULL) 22359b94acceSBarry Smith 22369b94acceSBarry Smith Notes: 22379b94acceSBarry Smith SNESGetFunction() is valid for SNES_NONLINEAR_EQUATIONS methods only 22389b94acceSBarry Smith Analogous routines for SNES_UNCONSTRAINED_MINIMIZATION methods are 22399b94acceSBarry Smith SNESGetMinimizationFunction() and SNESGetGradient(); 22409b94acceSBarry Smith 224136851e7fSLois Curfman McInnes Level: advanced 224236851e7fSLois Curfman McInnes 2243a86d99e1SLois Curfman McInnes .keywords: SNES, nonlinear, get, function 22449b94acceSBarry Smith 224531615d04SLois Curfman McInnes .seealso: SNESSetFunction(), SNESGetSolution(), SNESGetMinimizationFunction(), 224631615d04SLois Curfman McInnes SNESGetGradient() 22477bf4e008SBarry Smith 22489b94acceSBarry Smith @*/ 22497bf4e008SBarry Smith int SNESGetFunction(SNES snes,Vec *r,void **ctx) 22509b94acceSBarry Smith { 22513a40ed3dSBarry Smith PetscFunctionBegin; 225277c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2253a8c6a408SBarry Smith if (snes->method_class != SNES_NONLINEAR_EQUATIONS) { 2254a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_NONLINEAR_EQUATIONS only"); 2255a8c6a408SBarry Smith } 22567bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 22577bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 22583a40ed3dSBarry Smith PetscFunctionReturn(0); 22599b94acceSBarry Smith } 22609b94acceSBarry Smith 22615615d1e5SSatish Balay #undef __FUNC__ 2262d4bb536fSBarry Smith #define __FUNC__ "SNESGetGradient" 22639b94acceSBarry Smith /*@C 22643638b69dSLois Curfman McInnes SNESGetGradient - Returns the vector where the gradient is stored. 22659b94acceSBarry Smith 2266c7afd0dbSLois Curfman McInnes Not Collective, but Vec is parallel if SNES is parallel 2267c7afd0dbSLois Curfman McInnes 22689b94acceSBarry Smith Input Parameter: 22699b94acceSBarry Smith . snes - the SNES context 22709b94acceSBarry Smith 22719b94acceSBarry Smith Output Parameter: 22727bf4e008SBarry Smith + r - the gradient (or PETSC_NULL) 22737bf4e008SBarry Smith - ctx - the gradient context (or PETSC_NULL) 22749b94acceSBarry Smith 22759b94acceSBarry Smith Notes: 22769b94acceSBarry Smith SNESGetGradient() is valid for SNES_UNCONSTRAINED_MINIMIZATION methods 22779b94acceSBarry Smith only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 22789b94acceSBarry Smith SNESGetFunction(). 22799b94acceSBarry Smith 228036851e7fSLois Curfman McInnes Level: advanced 228136851e7fSLois Curfman McInnes 22829b94acceSBarry Smith .keywords: SNES, nonlinear, get, gradient 22839b94acceSBarry Smith 22847bf4e008SBarry Smith .seealso: SNESGetMinimizationFunction(), SNESGetSolution(), SNESGetFunction(), 22857bf4e008SBarry Smith SNESSetGradient(), SNESSetFunction() 22867bf4e008SBarry Smith 22879b94acceSBarry Smith @*/ 22887bf4e008SBarry Smith int SNESGetGradient(SNES snes,Vec *r,void **ctx) 22899b94acceSBarry Smith { 22903a40ed3dSBarry Smith PetscFunctionBegin; 229177c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 22923a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 2293a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 22943a40ed3dSBarry Smith } 22957bf4e008SBarry Smith if (r) *r = snes->vec_func_always; 22967bf4e008SBarry Smith if (ctx) *ctx = snes->funP; 22973a40ed3dSBarry Smith PetscFunctionReturn(0); 22989b94acceSBarry Smith } 22999b94acceSBarry Smith 23005615d1e5SSatish Balay #undef __FUNC__ 2301d4bb536fSBarry Smith #define __FUNC__ "SNESGetMinimizationFunction" 23027bf4e008SBarry Smith /*@C 23039b94acceSBarry Smith SNESGetMinimizationFunction - Returns the scalar function value for 23049b94acceSBarry Smith unconstrained minimization problems. 23059b94acceSBarry Smith 2306c7afd0dbSLois Curfman McInnes Not Collective 2307c7afd0dbSLois Curfman McInnes 23089b94acceSBarry Smith Input Parameter: 23099b94acceSBarry Smith . snes - the SNES context 23109b94acceSBarry Smith 23119b94acceSBarry Smith Output Parameter: 23127bf4e008SBarry Smith + r - the function (or PETSC_NULL) 23137bf4e008SBarry Smith - ctx - the function context (or PETSC_NULL) 23149b94acceSBarry Smith 23159b94acceSBarry Smith Notes: 23169b94acceSBarry Smith SNESGetMinimizationFunction() is valid for SNES_UNCONSTRAINED_MINIMIZATION 23179b94acceSBarry Smith methods only. An analogous routine for SNES_NONLINEAR_EQUATIONS methods is 23189b94acceSBarry Smith SNESGetFunction(). 23199b94acceSBarry Smith 232036851e7fSLois Curfman McInnes Level: advanced 232136851e7fSLois Curfman McInnes 23229b94acceSBarry Smith .keywords: SNES, nonlinear, get, function 23239b94acceSBarry Smith 23247bf4e008SBarry Smith .seealso: SNESGetGradient(), SNESGetSolution(), SNESGetFunction(), SNESSetFunction() 23257bf4e008SBarry Smith 23269b94acceSBarry Smith @*/ 23277bf4e008SBarry Smith int SNESGetMinimizationFunction(SNES snes,double *r,void **ctx) 23289b94acceSBarry Smith { 23293a40ed3dSBarry Smith PetscFunctionBegin; 233077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 233174679c65SBarry Smith PetscValidScalarPointer(r); 23323a40ed3dSBarry Smith if (snes->method_class != SNES_UNCONSTRAINED_MINIMIZATION) { 2333a8c6a408SBarry Smith SETERRQ(PETSC_ERR_ARG_WRONG,0,"For SNES_UNCONSTRAINED_MINIMIZATION only"); 23343a40ed3dSBarry Smith } 23357bf4e008SBarry Smith if (r) *r = snes->fc; 23367bf4e008SBarry Smith if (ctx) *ctx = snes->umfunP; 23373a40ed3dSBarry Smith PetscFunctionReturn(0); 23389b94acceSBarry Smith } 23399b94acceSBarry Smith 23405615d1e5SSatish Balay #undef __FUNC__ 2341d4bb536fSBarry Smith #define __FUNC__ "SNESSetOptionsPrefix" 23423c7409f5SSatish Balay /*@C 23433c7409f5SSatish Balay SNESSetOptionsPrefix - Sets the prefix used for searching for all 2344d850072dSLois Curfman McInnes SNES options in the database. 23453c7409f5SSatish Balay 2346fee21e36SBarry Smith Collective on SNES 2347fee21e36SBarry Smith 2348c7afd0dbSLois Curfman McInnes Input Parameter: 2349c7afd0dbSLois Curfman McInnes + snes - the SNES context 2350c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2351c7afd0dbSLois Curfman McInnes 2352d850072dSLois Curfman McInnes Notes: 2353a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2354c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2355d850072dSLois Curfman McInnes 235636851e7fSLois Curfman McInnes Level: advanced 235736851e7fSLois Curfman McInnes 23583c7409f5SSatish Balay .keywords: SNES, set, options, prefix, database 2359a86d99e1SLois Curfman McInnes 2360a86d99e1SLois Curfman McInnes .seealso: SNESSetFromOptions() 23613c7409f5SSatish Balay @*/ 23623c7409f5SSatish Balay int SNESSetOptionsPrefix(SNES snes,char *prefix) 23633c7409f5SSatish Balay { 23643c7409f5SSatish Balay int ierr; 23653c7409f5SSatish Balay 23663a40ed3dSBarry Smith PetscFunctionBegin; 236777c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2368639f9d9dSBarry Smith ierr = PetscObjectSetOptionsPrefix((PetscObject)snes, prefix);CHKERRQ(ierr); 23693c7409f5SSatish Balay ierr = SLESSetOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 23703a40ed3dSBarry Smith PetscFunctionReturn(0); 23713c7409f5SSatish Balay } 23723c7409f5SSatish Balay 23735615d1e5SSatish Balay #undef __FUNC__ 2374d4bb536fSBarry Smith #define __FUNC__ "SNESAppendOptionsPrefix" 23753c7409f5SSatish Balay /*@C 2376f525115eSLois Curfman McInnes SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2377d850072dSLois Curfman McInnes SNES options in the database. 23783c7409f5SSatish Balay 2379fee21e36SBarry Smith Collective on SNES 2380fee21e36SBarry Smith 2381c7afd0dbSLois Curfman McInnes Input Parameters: 2382c7afd0dbSLois Curfman McInnes + snes - the SNES context 2383c7afd0dbSLois Curfman McInnes - prefix - the prefix to prepend to all option names 2384c7afd0dbSLois Curfman McInnes 2385d850072dSLois Curfman McInnes Notes: 2386a83b1b31SSatish Balay A hyphen (-) must NOT be given at the beginning of the prefix name. 2387c7afd0dbSLois Curfman McInnes The first character of all runtime options is AUTOMATICALLY the hyphen. 2388d850072dSLois Curfman McInnes 238936851e7fSLois Curfman McInnes Level: advanced 239036851e7fSLois Curfman McInnes 23913c7409f5SSatish Balay .keywords: SNES, append, options, prefix, database 2392a86d99e1SLois Curfman McInnes 2393a86d99e1SLois Curfman McInnes .seealso: SNESGetOptionsPrefix() 23943c7409f5SSatish Balay @*/ 23953c7409f5SSatish Balay int SNESAppendOptionsPrefix(SNES snes,char *prefix) 23963c7409f5SSatish Balay { 23973c7409f5SSatish Balay int ierr; 23983c7409f5SSatish Balay 23993a40ed3dSBarry Smith PetscFunctionBegin; 240077c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2401639f9d9dSBarry Smith ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes, prefix);CHKERRQ(ierr); 24023c7409f5SSatish Balay ierr = SLESAppendOptionsPrefix(snes->sles,prefix);CHKERRQ(ierr); 24033a40ed3dSBarry Smith PetscFunctionReturn(0); 24043c7409f5SSatish Balay } 24053c7409f5SSatish Balay 24065615d1e5SSatish Balay #undef __FUNC__ 2407d4bb536fSBarry Smith #define __FUNC__ "SNESGetOptionsPrefix" 24089ab63eb5SSatish Balay /*@C 24093c7409f5SSatish Balay SNESGetOptionsPrefix - Sets the prefix used for searching for all 24103c7409f5SSatish Balay SNES options in the database. 24113c7409f5SSatish Balay 2412c7afd0dbSLois Curfman McInnes Not Collective 2413c7afd0dbSLois Curfman McInnes 24143c7409f5SSatish Balay Input Parameter: 24153c7409f5SSatish Balay . snes - the SNES context 24163c7409f5SSatish Balay 24173c7409f5SSatish Balay Output Parameter: 24183c7409f5SSatish Balay . prefix - pointer to the prefix string used 24193c7409f5SSatish Balay 24209ab63eb5SSatish Balay Notes: On the fortran side, the user should pass in a string 'prifix' of 24219ab63eb5SSatish Balay sufficient length to hold the prefix. 24229ab63eb5SSatish Balay 242336851e7fSLois Curfman McInnes Level: advanced 242436851e7fSLois Curfman McInnes 24253c7409f5SSatish Balay .keywords: SNES, get, options, prefix, database 2426a86d99e1SLois Curfman McInnes 2427a86d99e1SLois Curfman McInnes .seealso: SNESAppendOptionsPrefix() 24283c7409f5SSatish Balay @*/ 24293c7409f5SSatish Balay int SNESGetOptionsPrefix(SNES snes,char **prefix) 24303c7409f5SSatish Balay { 24313c7409f5SSatish Balay int ierr; 24323c7409f5SSatish Balay 24333a40ed3dSBarry Smith PetscFunctionBegin; 243477c4ece6SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2435639f9d9dSBarry Smith ierr = PetscObjectGetOptionsPrefix((PetscObject)snes, prefix);CHKERRQ(ierr); 24363a40ed3dSBarry Smith PetscFunctionReturn(0); 24373c7409f5SSatish Balay } 24383c7409f5SSatish Balay 2439ca161407SBarry Smith #undef __FUNC__ 2440ca161407SBarry Smith #define __FUNC__ "SNESPrintHelp" 2441ca161407SBarry Smith /*@ 2442ca161407SBarry Smith SNESPrintHelp - Prints all options for the SNES component. 2443ca161407SBarry Smith 2444c7afd0dbSLois Curfman McInnes Collective on SNES 2445c7afd0dbSLois Curfman McInnes 2446ca161407SBarry Smith Input Parameter: 2447ca161407SBarry Smith . snes - the SNES context 2448ca161407SBarry Smith 2449ca161407SBarry Smith Options Database Keys: 2450c7afd0dbSLois Curfman McInnes + -help - Prints SNES options 2451c7afd0dbSLois Curfman McInnes - -h - Prints SNES options 2452ca161407SBarry Smith 245336851e7fSLois Curfman McInnes Level: beginner 245436851e7fSLois Curfman McInnes 2455ca161407SBarry Smith .keywords: SNES, nonlinear, help 2456ca161407SBarry Smith 2457ca161407SBarry Smith .seealso: SNESSetFromOptions() 2458ca161407SBarry Smith @*/ 2459ca161407SBarry Smith int SNESPrintHelp(SNES snes) 2460ca161407SBarry Smith { 2461ca161407SBarry Smith char p[64]; 2462ca161407SBarry Smith SNES_KSP_EW_ConvCtx *kctx; 2463ca161407SBarry Smith int ierr; 2464ca161407SBarry Smith 2465ca161407SBarry Smith PetscFunctionBegin; 2466ca161407SBarry Smith PetscValidHeaderSpecific(snes,SNES_COOKIE); 2467ca161407SBarry Smith 2468549d3d68SSatish Balay ierr = PetscStrcpy(p,"-");CHKERRQ(ierr); 2469ca161407SBarry Smith if (snes->prefix) PetscStrcat(p, snes->prefix); 2470ca161407SBarry Smith 2471ca161407SBarry Smith kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 2472ca161407SBarry Smith 2473ebb8b11fSBarry Smith if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 2474d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm,"SNES options ------------------------------------------------\n");CHKERRQ(ierr); 2475488ecbafSBarry Smith ierr = FListPrintTypes(snes->comm,stdout,snes->prefix,"snes_type",SNESList);CHKERRQ(ierr); 2476d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_view: view SNES info after each nonlinear solve\n",p);CHKERRQ(ierr); 2477d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_max_it <its>: max iterations (default %d)\n",p,snes->max_its);CHKERRQ(ierr); 2478d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_max_funcs <maxf>: max function evals (default %d)\n",p,snes->max_funcs);CHKERRQ(ierr); 2479d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_stol <stol>: successive step tolerance (default %g)\n",p,snes->xtol);CHKERRQ(ierr); 2480d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_atol <atol>: absolute tolerance (default %g)\n",p,snes->atol);CHKERRQ(ierr); 2481d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_rtol <rtol>: relative tolerance (default %g)\n",p,snes->rtol);CHKERRQ(ierr); 2482d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_trtol <trtol>: trust region parameter tolerance (default %g)\n",p,snes->deltatol);CHKERRQ(ierr); 2483d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," SNES Monitoring Options: Choose any of the following\n");CHKERRQ(ierr); 2484d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_cancelmonitors: cancels all monitors hardwired in code\n",p);CHKERRQ(ierr); 2485d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_monitor: use default SNES convergence monitor, prints\n\ 2486d132466eSBarry Smith residual norm at each iteration.\n",p);CHKERRQ(ierr); 2487d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_smonitor: same as the above, but prints fewer digits of the\n\ 2488ca161407SBarry Smith residual norm for small residual norms. This is useful to conceal\n\ 2489d132466eSBarry Smith meaningless digits that may be different on different machines.\n",p);CHKERRQ(ierr); 2490d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_xmonitor [x,y,w,h]: use X graphics convergence monitor\n",p);CHKERRQ(ierr); 2491d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_vecmonitor: plots solution at each iteration \n",p);CHKERRQ(ierr); 2492d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_vecmonitor_update: plots update to solution at each iteration \n",p);CHKERRQ(ierr); 2493ca161407SBarry Smith if (snes->type == SNES_NONLINEAR_EQUATIONS) { 2494d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2495d132466eSBarry Smith " Options for solving systems of nonlinear equations only:\n");CHKERRQ(ierr); 2496d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Jacobian\n",p);CHKERRQ(ierr); 2497d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Jacobian\n",p);CHKERRQ(ierr); 2498d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_mf_operator:use matrix-free Jacobian and user-provided preconditioning matrix\n",p);CHKERRQ(ierr); 2499d132466eSBarry 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); 2500d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_no_convergence_test: Do not test for convergence, always run to SNES max its\n",p);CHKERRQ(ierr); 2501d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_ksp_ew_conv: use Eisenstat-Walker computation of KSP rtol. Params are:\n",p);CHKERRQ(ierr); 2502d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2503d132466eSBarry Smith " %ssnes_ksp_ew_version <version> (1 or 2, default is %d)\n",p,kctx->version);CHKERRQ(ierr); 2504d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2505d132466eSBarry Smith " %ssnes_ksp_ew_rtol0 <rtol0> (0 <= rtol0 < 1, default %g)\n",p,kctx->rtol_0);CHKERRQ(ierr); 2506d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2507d132466eSBarry Smith " %ssnes_ksp_ew_rtolmax <rtolmax> (0 <= rtolmax < 1, default %g)\n",p,kctx->rtol_max);CHKERRQ(ierr); 2508d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2509d132466eSBarry Smith " %ssnes_ksp_ew_gamma <gamma> (0 <= gamma <= 1, default %g)\n",p,kctx->gamma);CHKERRQ(ierr); 2510d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2511d132466eSBarry Smith " %ssnes_ksp_ew_alpha <alpha> (1 < alpha <= 2, default %g)\n",p,kctx->alpha);CHKERRQ(ierr); 2512d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2513d132466eSBarry Smith " %ssnes_ksp_ew_alpha2 <alpha2> (default %g)\n",p,kctx->alpha2);CHKERRQ(ierr); 2514d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm, 2515d132466eSBarry Smith " %ssnes_ksp_ew_threshold <threshold> (0 < threshold < 1, default %g)\n",p,kctx->threshold);CHKERRQ(ierr); 2516ca161407SBarry Smith } else if (snes->type == SNES_UNCONSTRAINED_MINIMIZATION) { 2517d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," Options for solving unconstrained minimization problems only:\n");CHKERRQ(ierr); 2518d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_fmin <ftol>: minimum function tolerance (default %g)\n",p,snes->fmin);CHKERRQ(ierr); 2519d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_fd: use finite differences for Hessian\n",p);CHKERRQ(ierr); 2520d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," %ssnes_mf: use matrix-free Hessian\n",p);CHKERRQ(ierr); 2521ca161407SBarry Smith } 2522454a90a3SBarry Smith ierr = (*PetscHelpPrintf)(snes->comm," Run program with -help %ssnes_type <type> for help on ",p);CHKERRQ(ierr); 2523d132466eSBarry Smith ierr = (*PetscHelpPrintf)(snes->comm,"a particular method\n");CHKERRQ(ierr); 2524ca161407SBarry Smith if (snes->printhelp) { 2525ca161407SBarry Smith ierr = (*snes->printhelp)(snes,p);CHKERRQ(ierr); 2526ca161407SBarry Smith } 2527ca161407SBarry Smith PetscFunctionReturn(0); 2528ca161407SBarry Smith } 25293c7409f5SSatish Balay 2530acb85ae6SSatish Balay /*MC 2531f1af5d2fSBarry Smith SNESRegisterDynamic - Adds a method to the nonlinear solver package. 25329b94acceSBarry Smith 2533b2002411SLois Curfman McInnes Synopsis: 2534f1af5d2fSBarry Smith SNESRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(SNES)) 25359b94acceSBarry Smith 25368d76a1e5SLois Curfman McInnes Not collective 25378d76a1e5SLois Curfman McInnes 2538b2002411SLois Curfman McInnes Input Parameters: 2539c7afd0dbSLois Curfman McInnes + name_solver - name of a new user-defined solver 2540b2002411SLois Curfman McInnes . path - path (either absolute or relative) the library containing this solver 2541b2002411SLois Curfman McInnes . name_create - name of routine to create method context 2542c7afd0dbSLois Curfman McInnes - routine_create - routine to create method context 25439b94acceSBarry Smith 2544b2002411SLois Curfman McInnes Notes: 2545f1af5d2fSBarry Smith SNESRegisterDynamic() may be called multiple times to add several user-defined solvers. 25463a40ed3dSBarry Smith 2547b2002411SLois Curfman McInnes If dynamic libraries are used, then the fourth input argument (routine_create) 2548b2002411SLois Curfman McInnes is ignored. 2549b2002411SLois Curfman McInnes 2550b2002411SLois Curfman McInnes Sample usage: 255169225d78SLois Curfman McInnes .vb 2552f1af5d2fSBarry Smith SNESRegisterDynamic("my_solver",/home/username/my_lib/lib/libg/solaris/mylib.a, 2553b2002411SLois Curfman McInnes "MySolverCreate",MySolverCreate); 255469225d78SLois Curfman McInnes .ve 2555b2002411SLois Curfman McInnes 2556b2002411SLois Curfman McInnes Then, your solver can be chosen with the procedural interface via 2557b2002411SLois Curfman McInnes $ SNESSetType(snes,"my_solver") 2558b2002411SLois Curfman McInnes or at runtime via the option 2559b2002411SLois Curfman McInnes $ -snes_type my_solver 2560b2002411SLois Curfman McInnes 256136851e7fSLois Curfman McInnes Level: advanced 256236851e7fSLois Curfman McInnes 256385614651SBarry Smith $PETSC_ARCH, $PETSC_DIR, $PETSC_LDIR, and $BOPT occuring in pathname will be replaced with appropriate values. 2564dd438238SBarry Smith 2565b2002411SLois Curfman McInnes .keywords: SNES, nonlinear, register 2566b2002411SLois Curfman McInnes 2567b2002411SLois Curfman McInnes .seealso: SNESRegisterAll(), SNESRegisterDestroy() 2568b2002411SLois Curfman McInnes M*/ 2569b2002411SLois Curfman McInnes 2570b2002411SLois Curfman McInnes #undef __FUNC__ 2571f1af5d2fSBarry Smith #define __FUNC__ "SNESRegister" 2572f1af5d2fSBarry Smith int SNESRegister(char *sname,char *path,char *name,int (*function)(SNES)) 2573b2002411SLois Curfman McInnes { 2574b2002411SLois Curfman McInnes char fullname[256]; 2575b2002411SLois Curfman McInnes int ierr; 2576b2002411SLois Curfman McInnes 2577b2002411SLois Curfman McInnes PetscFunctionBegin; 2578f1af5d2fSBarry Smith ierr = FListConcat(path,name,fullname); CHKERRQ(ierr); 2579f1af5d2fSBarry Smith ierr = FListAdd(&SNESList,sname,fullname, (int (*)(void*))function);CHKERRQ(ierr); 2580b2002411SLois Curfman McInnes PetscFunctionReturn(0); 2581b2002411SLois Curfman McInnes } 2582