1 #define PETSCSNES_DLL 2 3 #include "include/private/snesimpl.h" /*I "petscsnes.h" I*/ 4 5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 6 PetscFList SNESList = PETSC_NULL; 7 8 /* Logging support */ 9 PetscCookie PETSCSNES_DLLEXPORT SNES_COOKIE; 10 PetscLogEvent SNES_Solve, SNES_LineSearch, SNES_FunctionEval, SNES_JacobianEval; 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "SNESSetFunctionDomainError" 14 /*@ 15 SNESSetFunctionDomainError - tells SNES that the input vector to your FormFunction is not 16 in the functions domain. For example, negative pressure. 17 18 Collective on SNES 19 20 Input Parameters: 21 . SNES - the SNES context 22 23 Level: advanced 24 25 .keywords: SNES, view 26 27 .seealso: SNESCreate(), SNESSetFunction() 28 @*/ 29 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunctionDomainError(SNES snes) 30 { 31 PetscFunctionBegin; 32 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 33 snes->domainerror = PETSC_TRUE; 34 PetscFunctionReturn(0); 35 } 36 37 #undef __FUNCT__ 38 #define __FUNCT__ "SNESView" 39 /*@C 40 SNESView - Prints the SNES data structure. 41 42 Collective on SNES 43 44 Input Parameters: 45 + SNES - the SNES context 46 - viewer - visualization context 47 48 Options Database Key: 49 . -snes_view - Calls SNESView() at end of SNESSolve() 50 51 Notes: 52 The available visualization contexts include 53 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 54 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 55 output where only the first processor opens 56 the file. All other processors send their 57 data to the first processor to print. 58 59 The user can open an alternative visualization context with 60 PetscViewerASCIIOpen() - output to a specified file. 61 62 Level: beginner 63 64 .keywords: SNES, view 65 66 .seealso: PetscViewerASCIIOpen() 67 @*/ 68 PetscErrorCode PETSCSNES_DLLEXPORT SNESView(SNES snes,PetscViewer viewer) 69 { 70 SNESKSPEW *kctx; 71 PetscErrorCode ierr; 72 KSP ksp; 73 const SNESType type; 74 PetscTruth iascii,isstring; 75 76 PetscFunctionBegin; 77 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 78 if (!viewer) { 79 ierr = PetscViewerASCIIGetStdout(((PetscObject)snes)->comm,&viewer);CHKERRQ(ierr); 80 } 81 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 82 PetscCheckSameComm(snes,1,viewer,2); 83 84 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 85 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 86 if (iascii) { 87 if (((PetscObject)snes)->prefix) { 88 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",((PetscObject)snes)->prefix);CHKERRQ(ierr); 89 } else { 90 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 91 } 92 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 93 if (type) { 94 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 95 } else { 96 ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 97 } 98 if (snes->ops->view) { 99 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 100 ierr = (*snes->ops->view)(snes,viewer);CHKERRQ(ierr); 101 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 102 } 103 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 104 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%G, absolute=%G, solution=%G\n", 105 snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr); 106 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 107 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 108 if (snes->ksp_ewconv) { 109 kctx = (SNESKSPEW *)snes->kspconvctx; 110 if (kctx) { 111 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 112 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%G, rtol_max=%G, threshold=%G\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 113 ierr = PetscViewerASCIIPrintf(viewer," gamma=%G, alpha=%G, alpha2=%G\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 114 } 115 } 116 } else if (isstring) { 117 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 118 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 119 } 120 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 121 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 122 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 123 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 124 PetscFunctionReturn(0); 125 } 126 127 /* 128 We retain a list of functions that also take SNES command 129 line options. These are called at the end SNESSetFromOptions() 130 */ 131 #define MAXSETFROMOPTIONS 5 132 static PetscInt numberofsetfromoptions; 133 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 134 135 #undef __FUNCT__ 136 #define __FUNCT__ "SNESAddOptionsChecker" 137 /*@C 138 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 139 140 Not Collective 141 142 Input Parameter: 143 . snescheck - function that checks for options 144 145 Level: developer 146 147 .seealso: SNESSetFromOptions() 148 @*/ 149 PetscErrorCode PETSCSNES_DLLEXPORT SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 150 { 151 PetscFunctionBegin; 152 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 153 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 154 } 155 othersetfromoptions[numberofsetfromoptions++] = snescheck; 156 PetscFunctionReturn(0); 157 } 158 159 #undef __FUNCT__ 160 #define __FUNCT__ "SNESSetFromOptions" 161 /*@ 162 SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 163 164 Collective on SNES 165 166 Input Parameter: 167 . snes - the SNES context 168 169 Options Database Keys: 170 + -snes_type <type> - ls, tr, umls, umtr, test 171 . -snes_stol - convergence tolerance in terms of the norm 172 of the change in the solution between steps 173 . -snes_atol <abstol> - absolute tolerance of residual norm 174 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 175 . -snes_max_it <max_it> - maximum number of iterations 176 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 177 . -snes_max_fail <max_fail> - maximum number of failures 178 . -snes_max_linear_solve_fail - number of linear solver failures before SNESSolve() stops 179 . -snes_lag_preconditioner <lag> - how often preconditioner is rebuilt (use -1 to never rebuild) 180 . -snes_lag_jacobian <lag> - how often Jacobian is rebuilt (use -1 to never rebuild) 181 . -snes_trtol <trtol> - trust region tolerance 182 . -snes_no_convergence_test - skip convergence test in nonlinear 183 solver; hence iterations will continue until max_it 184 or some other criterion is reached. Saves expense 185 of convergence test 186 . -snes_monitor <optional filename> - prints residual norm at each iteration. if no 187 filename given prints to stdout 188 . -snes_monitor_solution - plots solution at each iteration 189 . -snes_monitor_residual - plots residual (not its norm) at each iteration 190 . -snes_monitor_solution_update - plots update to solution at each iteration 191 . -snes_monitor_draw - plots residual norm at each iteration 192 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 193 . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 194 - -snes_converged_reason - print the reason for convergence/divergence after each solve 195 196 Options Database for Eisenstat-Walker method: 197 + -snes_ksp_ew - use Eisenstat-Walker method for determining linear system convergence 198 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 199 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 200 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 201 . -snes_ksp_ew_gamma <gamma> - Sets gamma 202 . -snes_ksp_ew_alpha <alpha> - Sets alpha 203 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 204 - -snes_ksp_ew_threshold <threshold> - Sets threshold 205 206 Notes: 207 To see all options, run your program with the -help option or consult 208 the users manual. 209 210 Level: beginner 211 212 .keywords: SNES, nonlinear, set, options, database 213 214 .seealso: SNESSetOptionsPrefix() 215 @*/ 216 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFromOptions(SNES snes) 217 { 218 PetscTruth flg; 219 PetscInt i,indx,lag; 220 const char *deft = SNESLS; 221 const char *convtests[] = {"default","skip"}; 222 SNESKSPEW *kctx = NULL; 223 char type[256], monfilename[PETSC_MAX_PATH_LEN]; 224 PetscViewerASCIIMonitor monviewer; 225 PetscErrorCode ierr; 226 227 PetscFunctionBegin; 228 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 229 230 ierr = PetscOptionsBegin(((PetscObject)snes)->comm,((PetscObject)snes)->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 231 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 232 if (((PetscObject)snes)->type_name) { deft = ((PetscObject)snes)->type_name; } 233 ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 234 if (flg) { 235 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 236 } else if (!((PetscObject)snes)->type_name) { 237 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 238 } 239 ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 240 241 ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 242 ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr); 243 244 ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 245 ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 246 ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 247 ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 248 ierr = PetscOptionsInt("-snes_max_linear_solve_fail","Maximum failures in linear solves allowed","SNESSetMaxLinearSolveFailures",snes->maxLinearSolveFailures,&snes->maxLinearSolveFailures,PETSC_NULL);CHKERRQ(ierr); 249 250 ierr = PetscOptionsInt("-snes_lag_preconditioner","How often to rebuild preconditioner","SNESSetLagPreconditioner",snes->lagpreconditioner,&lag,&flg);CHKERRQ(ierr); 251 if (flg) { 252 ierr = SNESSetLagPreconditioner(snes,lag);CHKERRQ(ierr); 253 } 254 ierr = PetscOptionsInt("-snes_lag_jacobian","How often to rebuild Jacobian","SNESSetLagJacobian",snes->lagjacobian,&lag,&flg);CHKERRQ(ierr); 255 if (flg) { 256 ierr = SNESSetLagJacobian(snes,lag);CHKERRQ(ierr); 257 } 258 259 ierr = PetscOptionsEList("-snes_convergence_test","Convergence test","SNESSetConvergenceTest",convtests,2,"default",&indx,&flg);CHKERRQ(ierr); 260 if (flg) { 261 switch (indx) { 262 case 0: ierr = SNESSetConvergenceTest(snes,SNESDefaultConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 263 case 1: ierr = SNESSetConvergenceTest(snes,SNESSkipConverged,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr); break; 264 } 265 } 266 267 ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr); 268 if (flg) { 269 snes->printreason = PETSC_TRUE; 270 } 271 272 kctx = (SNESKSPEW *)snes->kspconvctx; 273 274 ierr = PetscOptionsTruth("-snes_ksp_ew","Use Eisentat-Walker linear system convergence test","SNESKSPSetUseEW",snes->ksp_ewconv,&snes->ksp_ewconv,PETSC_NULL);CHKERRQ(ierr); 275 276 ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1, 2 or 3","SNESKSPSetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 277 ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNESKSPSetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 278 ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNESKSPSetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 279 ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNESKSPSetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 280 ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNESKSPSetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 281 ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNESKSPSetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 282 ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNESKSPSetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 283 284 ierr = PetscOptionsName("-snes_monitor_cancel","Remove all monitors","SNESMonitorCancel",&flg);CHKERRQ(ierr); 285 if (flg) {ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);} 286 287 ierr = PetscOptionsString("-snes_monitor","Monitor norm of function","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 288 if (flg) { 289 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 290 ierr = SNESMonitorSet(snes,SNESMonitorDefault,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 291 } 292 293 ierr = PetscOptionsString("-snes_ratiomonitor","Monitor ratios of norms of function","SNESMonitorSetRatio","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 294 if (flg) { 295 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 296 ierr = SNESMonitorSetRatio(snes,monviewer);CHKERRQ(ierr); 297 } 298 299 ierr = PetscOptionsString("-snes_monitor_short","Monitor norm of function (fewer digits)","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 300 if (flg) { 301 ierr = PetscViewerASCIIMonitorCreate(((PetscObject)snes)->comm,monfilename,((PetscObject)snes)->tablevel,&monviewer);CHKERRQ(ierr); 302 ierr = SNESMonitorSet(snes,SNESMonitorDefaultShort,monviewer,(PetscErrorCode (*)(void*))PetscViewerASCIIMonitorDestroy);CHKERRQ(ierr); 303 } 304 305 ierr = PetscOptionsName("-snes_monitor_solution","Plot solution at each iteration","SNESMonitorSolution",&flg);CHKERRQ(ierr); 306 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolution,0,0);CHKERRQ(ierr);} 307 ierr = PetscOptionsName("-snes_monitor_solution_update","Plot correction at each iteration","SNESMonitorSolutionUpdate",&flg);CHKERRQ(ierr); 308 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorSolutionUpdate,0,0);CHKERRQ(ierr);} 309 ierr = PetscOptionsName("-snes_monitor_residual","Plot residual at each iteration","SNESMonitorResidual",&flg);CHKERRQ(ierr); 310 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorResidual,0,0);CHKERRQ(ierr);} 311 ierr = PetscOptionsName("-snes_monitor_draw","Plot function norm at each iteration","SNESMonitorLG",&flg);CHKERRQ(ierr); 312 if (flg) {ierr = SNESMonitorSet(snes,SNESMonitorLG,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 313 314 ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 315 if (flg) { 316 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 317 ierr = PetscInfo(snes,"Setting default finite difference Jacobian matrix\n");CHKERRQ(ierr); 318 } 319 320 for(i = 0; i < numberofsetfromoptions; i++) { 321 ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 322 } 323 324 if (snes->ops->setfromoptions) { 325 ierr = (*snes->ops->setfromoptions)(snes);CHKERRQ(ierr); 326 } 327 ierr = PetscOptionsEnd();CHKERRQ(ierr); 328 329 330 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 331 ierr = KSPSetFromOptions(snes->ksp);CHKERRQ(ierr); 332 333 PetscFunctionReturn(0); 334 } 335 336 337 #undef __FUNCT__ 338 #define __FUNCT__ "SNESSetApplicationContext" 339 /*@ 340 SNESSetApplicationContext - Sets the optional user-defined context for 341 the nonlinear solvers. 342 343 Collective on SNES 344 345 Input Parameters: 346 + snes - the SNES context 347 - usrP - optional user context 348 349 Level: intermediate 350 351 .keywords: SNES, nonlinear, set, application, context 352 353 .seealso: SNESGetApplicationContext() 354 @*/ 355 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetApplicationContext(SNES snes,void *usrP) 356 { 357 PetscFunctionBegin; 358 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 359 snes->user = usrP; 360 PetscFunctionReturn(0); 361 } 362 363 #undef __FUNCT__ 364 #define __FUNCT__ "SNESGetApplicationContext" 365 /*@C 366 SNESGetApplicationContext - Gets the user-defined context for the 367 nonlinear solvers. 368 369 Not Collective 370 371 Input Parameter: 372 . snes - SNES context 373 374 Output Parameter: 375 . usrP - user context 376 377 Level: intermediate 378 379 .keywords: SNES, nonlinear, get, application, context 380 381 .seealso: SNESSetApplicationContext() 382 @*/ 383 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetApplicationContext(SNES snes,void **usrP) 384 { 385 PetscFunctionBegin; 386 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 387 *usrP = snes->user; 388 PetscFunctionReturn(0); 389 } 390 391 #undef __FUNCT__ 392 #define __FUNCT__ "SNESGetIterationNumber" 393 /*@ 394 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 395 at this time. 396 397 Not Collective 398 399 Input Parameter: 400 . snes - SNES context 401 402 Output Parameter: 403 . iter - iteration number 404 405 Notes: 406 For example, during the computation of iteration 2 this would return 1. 407 408 This is useful for using lagged Jacobians (where one does not recompute the 409 Jacobian at each SNES iteration). For example, the code 410 .vb 411 ierr = SNESGetIterationNumber(snes,&it); 412 if (!(it % 2)) { 413 [compute Jacobian here] 414 } 415 .ve 416 can be used in your ComputeJacobian() function to cause the Jacobian to be 417 recomputed every second SNES iteration. 418 419 Level: intermediate 420 421 .keywords: SNES, nonlinear, get, iteration, number, 422 423 .seealso: SNESGetFunctionNorm(), SNESGetLinearSolveIterations() 424 @*/ 425 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetIterationNumber(SNES snes,PetscInt* iter) 426 { 427 PetscFunctionBegin; 428 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 429 PetscValidIntPointer(iter,2); 430 *iter = snes->iter; 431 PetscFunctionReturn(0); 432 } 433 434 #undef __FUNCT__ 435 #define __FUNCT__ "SNESGetFunctionNorm" 436 /*@ 437 SNESGetFunctionNorm - Gets the norm of the current function that was set 438 with SNESSSetFunction(). 439 440 Collective on SNES 441 442 Input Parameter: 443 . snes - SNES context 444 445 Output Parameter: 446 . fnorm - 2-norm of function 447 448 Level: intermediate 449 450 .keywords: SNES, nonlinear, get, function, norm 451 452 .seealso: SNESGetFunction(), SNESGetIterationNumber(), SNESGetLinearSolveIterations() 453 @*/ 454 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunctionNorm(SNES snes,PetscReal *fnorm) 455 { 456 PetscFunctionBegin; 457 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 458 PetscValidScalarPointer(fnorm,2); 459 *fnorm = snes->norm; 460 PetscFunctionReturn(0); 461 } 462 463 #undef __FUNCT__ 464 #define __FUNCT__ "SNESGetNonlinearStepFailures" 465 /*@ 466 SNESGetNonlinearStepFailures - Gets the number of unsuccessful steps 467 attempted by the nonlinear solver. 468 469 Not Collective 470 471 Input Parameter: 472 . snes - SNES context 473 474 Output Parameter: 475 . nfails - number of unsuccessful steps attempted 476 477 Notes: 478 This counter is reset to zero for each successive call to SNESSolve(). 479 480 Level: intermediate 481 482 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 483 484 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 485 SNESSetMaxNonlinearStepFailures(), SNESGetMaxNonlinearStepFailures() 486 @*/ 487 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNonlinearStepFailures(SNES snes,PetscInt* nfails) 488 { 489 PetscFunctionBegin; 490 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 491 PetscValidIntPointer(nfails,2); 492 *nfails = snes->numFailures; 493 PetscFunctionReturn(0); 494 } 495 496 #undef __FUNCT__ 497 #define __FUNCT__ "SNESSetMaxNonlinearStepFailures" 498 /*@ 499 SNESSetMaxNonlinearStepFailures - Sets the maximum number of unsuccessful steps 500 attempted by the nonlinear solver before it gives up. 501 502 Not Collective 503 504 Input Parameters: 505 + snes - SNES context 506 - maxFails - maximum of unsuccessful steps 507 508 Level: intermediate 509 510 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 511 512 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 513 SNESGetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 514 @*/ 515 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxNonlinearStepFailures(SNES snes, PetscInt maxFails) 516 { 517 PetscFunctionBegin; 518 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 519 snes->maxFailures = maxFails; 520 PetscFunctionReturn(0); 521 } 522 523 #undef __FUNCT__ 524 #define __FUNCT__ "SNESGetMaxNonlinearStepFailures" 525 /*@ 526 SNESGetMaxNonlinearStepFailures - Gets the maximum number of unsuccessful steps 527 attempted by the nonlinear solver before it gives up. 528 529 Not Collective 530 531 Input Parameter: 532 . snes - SNES context 533 534 Output Parameter: 535 . maxFails - maximum of unsuccessful steps 536 537 Level: intermediate 538 539 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 540 541 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures(), 542 SNESSetMaxNonlinearStepFailures(), SNESGetNonlinearStepFailures() 543 544 @*/ 545 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxNonlinearStepFailures(SNES snes, PetscInt *maxFails) 546 { 547 PetscFunctionBegin; 548 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 549 PetscValidIntPointer(maxFails,2); 550 *maxFails = snes->maxFailures; 551 PetscFunctionReturn(0); 552 } 553 554 #undef __FUNCT__ 555 #define __FUNCT__ "SNESGetNumberFunctionEvals" 556 /*@ 557 SNESGetNumberFunctionEvals - Gets the number of user provided function evaluations 558 done by SNES. 559 560 Not Collective 561 562 Input Parameter: 563 . snes - SNES context 564 565 Output Parameter: 566 . nfuncs - number of evaluations 567 568 Level: intermediate 569 570 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 571 572 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), SNESGetLinearSolveFailures() 573 @*/ 574 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetNumberFunctionEvals(SNES snes, PetscInt *nfuncs) 575 { 576 PetscFunctionBegin; 577 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 578 PetscValidIntPointer(nfuncs,2); 579 *nfuncs = snes->nfuncs; 580 PetscFunctionReturn(0); 581 } 582 583 #undef __FUNCT__ 584 #define __FUNCT__ "SNESGetLinearSolveFailures" 585 /*@ 586 SNESGetLinearSolveFailures - Gets the number of failed (non-converged) 587 linear solvers. 588 589 Not Collective 590 591 Input Parameter: 592 . snes - SNES context 593 594 Output Parameter: 595 . nfails - number of failed solves 596 597 Notes: 598 This counter is reset to zero for each successive call to SNESSolve(). 599 600 Level: intermediate 601 602 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 603 604 .seealso: SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures() 605 @*/ 606 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveFailures(SNES snes,PetscInt* nfails) 607 { 608 PetscFunctionBegin; 609 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 610 PetscValidIntPointer(nfails,2); 611 *nfails = snes->numLinearSolveFailures; 612 PetscFunctionReturn(0); 613 } 614 615 #undef __FUNCT__ 616 #define __FUNCT__ "SNESSetMaxLinearSolveFailures" 617 /*@ 618 SNESSetMaxLinearSolveFailures - the number of failed linear solve attempts 619 allowed before SNES returns with a diverged reason of SNES_DIVERGED_LINEAR_SOLVE 620 621 Collective on SNES 622 623 Input Parameters: 624 + snes - SNES context 625 - maxFails - maximum allowed linear solve failures 626 627 Level: intermediate 628 629 Notes: By default this is 0; that is SNES returns on the first failed linear solve 630 631 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 632 633 .seealso: SNESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures(), SNESGetLinearSolveIterations() 634 @*/ 635 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetMaxLinearSolveFailures(SNES snes, PetscInt maxFails) 636 { 637 PetscFunctionBegin; 638 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 639 snes->maxLinearSolveFailures = maxFails; 640 PetscFunctionReturn(0); 641 } 642 643 #undef __FUNCT__ 644 #define __FUNCT__ "SNESGetMaxLinearSolveFailures" 645 /*@ 646 SNESGetMaxLinearSolveFailures - gets the maximum number of linear solve failures that 647 are allowed before SNES terminates 648 649 Not Collective 650 651 Input Parameter: 652 . snes - SNES context 653 654 Output Parameter: 655 . maxFails - maximum of unsuccessful solves allowed 656 657 Level: intermediate 658 659 Notes: By default this is 1; that is SNES returns on the first failed linear solve 660 661 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 662 663 .seealso: SNESGetLinearSolveFailures(), SNESGetLinearSolveIterations(), SNESSetMaxLinearSolveFailures(), 664 @*/ 665 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetMaxLinearSolveFailures(SNES snes, PetscInt *maxFails) 666 { 667 PetscFunctionBegin; 668 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 669 PetscValidIntPointer(maxFails,2); 670 *maxFails = snes->maxLinearSolveFailures; 671 PetscFunctionReturn(0); 672 } 673 674 #undef __FUNCT__ 675 #define __FUNCT__ "SNESGetLinearSolveIterations" 676 /*@ 677 SNESGetLinearSolveIterations - Gets the total number of linear iterations 678 used by the nonlinear solver. 679 680 Not Collective 681 682 Input Parameter: 683 . snes - SNES context 684 685 Output Parameter: 686 . lits - number of linear iterations 687 688 Notes: 689 This counter is reset to zero for each successive call to SNESSolve(). 690 691 Level: intermediate 692 693 .keywords: SNES, nonlinear, get, number, linear, iterations 694 695 .seealso: SNESGetIterationNumber(), SNESGetFunctionNorm()S, NESGetLinearSolveFailures(), SNESGetMaxLinearSolveFailures() 696 @*/ 697 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLinearSolveIterations(SNES snes,PetscInt* lits) 698 { 699 PetscFunctionBegin; 700 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 701 PetscValidIntPointer(lits,2); 702 *lits = snes->linear_its; 703 PetscFunctionReturn(0); 704 } 705 706 #undef __FUNCT__ 707 #define __FUNCT__ "SNESGetKSP" 708 /*@ 709 SNESGetKSP - Returns the KSP context for a SNES solver. 710 711 Not Collective, but if SNES object is parallel, then KSP object is parallel 712 713 Input Parameter: 714 . snes - the SNES context 715 716 Output Parameter: 717 . ksp - the KSP context 718 719 Notes: 720 The user can then directly manipulate the KSP context to set various 721 options, etc. Likewise, the user can then extract and manipulate the 722 PC contexts as well. 723 724 Level: beginner 725 726 .keywords: SNES, nonlinear, get, KSP, context 727 728 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 729 @*/ 730 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetKSP(SNES snes,KSP *ksp) 731 { 732 PetscErrorCode ierr; 733 734 PetscFunctionBegin; 735 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 736 PetscValidPointer(ksp,2); 737 738 if (!snes->ksp) { 739 ierr = KSPCreate(((PetscObject)snes)->comm,&snes->ksp);CHKERRQ(ierr); 740 ierr = PetscObjectIncrementTabLevel((PetscObject)snes->ksp,(PetscObject)snes,1);CHKERRQ(ierr); 741 ierr = PetscLogObjectParent(snes,snes->ksp);CHKERRQ(ierr); 742 } 743 *ksp = snes->ksp; 744 PetscFunctionReturn(0); 745 } 746 747 #undef __FUNCT__ 748 #define __FUNCT__ "SNESSetKSP" 749 /*@ 750 SNESSetKSP - Sets a KSP context for the SNES object to use 751 752 Not Collective, but the SNES and KSP objects must live on the same MPI_Comm 753 754 Input Parameters: 755 + snes - the SNES context 756 - ksp - the KSP context 757 758 Notes: 759 The SNES object already has its KSP object, you can obtain with SNESGetKSP() 760 so this routine is rarely needed. 761 762 The KSP object that is already in the SNES object has its reference count 763 decreased by one. 764 765 Level: developer 766 767 .keywords: SNES, nonlinear, get, KSP, context 768 769 .seealso: KSPGetPC(), SNESCreate(), KSPCreate(), SNESSetKSP() 770 @*/ 771 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetKSP(SNES snes,KSP ksp) 772 { 773 PetscErrorCode ierr; 774 775 PetscFunctionBegin; 776 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 777 PetscValidHeaderSpecific(ksp,KSP_COOKIE,2); 778 PetscCheckSameComm(snes,1,ksp,2); 779 ierr = PetscObjectReference((PetscObject)ksp);CHKERRQ(ierr); 780 if (snes->ksp) {ierr = PetscObjectDereference((PetscObject)snes->ksp);CHKERRQ(ierr);} 781 snes->ksp = ksp; 782 PetscFunctionReturn(0); 783 } 784 785 #if 0 786 #undef __FUNCT__ 787 #define __FUNCT__ "SNESPublish_Petsc" 788 static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 789 { 790 PetscFunctionBegin; 791 PetscFunctionReturn(0); 792 } 793 #endif 794 795 /* -----------------------------------------------------------*/ 796 #undef __FUNCT__ 797 #define __FUNCT__ "SNESCreate" 798 /*@ 799 SNESCreate - Creates a nonlinear solver context. 800 801 Collective on MPI_Comm 802 803 Input Parameters: 804 . comm - MPI communicator 805 806 Output Parameter: 807 . outsnes - the new SNES context 808 809 Options Database Keys: 810 + -snes_mf - Activates default matrix-free Jacobian-vector products, 811 and no preconditioning matrix 812 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 813 products, and a user-provided preconditioning matrix 814 as set by SNESSetJacobian() 815 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 816 817 Level: beginner 818 819 .keywords: SNES, nonlinear, create, context 820 821 .seealso: SNESSolve(), SNESDestroy(), SNES, SNESSetLagPreconditioner() 822 823 @*/ 824 PetscErrorCode PETSCSNES_DLLEXPORT SNESCreate(MPI_Comm comm,SNES *outsnes) 825 { 826 PetscErrorCode ierr; 827 SNES snes; 828 SNESKSPEW *kctx; 829 830 PetscFunctionBegin; 831 PetscValidPointer(outsnes,2); 832 *outsnes = PETSC_NULL; 833 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 834 ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 835 #endif 836 837 ierr = PetscHeaderCreate(snes,_p_SNES,struct _SNESOps,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView);CHKERRQ(ierr); 838 839 snes->ops->converged = SNESDefaultConverged; 840 snes->max_its = 50; 841 snes->max_funcs = 10000; 842 snes->norm = 0.0; 843 snes->rtol = 1.e-8; 844 snes->ttol = 0.0; 845 snes->abstol = 1.e-50; 846 snes->xtol = 1.e-8; 847 snes->deltatol = 1.e-12; 848 snes->nfuncs = 0; 849 snes->numFailures = 0; 850 snes->maxFailures = 1; 851 snes->linear_its = 0; 852 snes->lagjacobian = 1; 853 snes->lagpreconditioner = 1; 854 snes->numbermonitors = 0; 855 snes->data = 0; 856 snes->setupcalled = PETSC_FALSE; 857 snes->ksp_ewconv = PETSC_FALSE; 858 snes->vwork = 0; 859 snes->nwork = 0; 860 snes->conv_hist_len = 0; 861 snes->conv_hist_max = 0; 862 snes->conv_hist = PETSC_NULL; 863 snes->conv_hist_its = PETSC_NULL; 864 snes->conv_hist_reset = PETSC_TRUE; 865 snes->reason = SNES_CONVERGED_ITERATING; 866 867 snes->numLinearSolveFailures = 0; 868 snes->maxLinearSolveFailures = 1; 869 870 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 871 ierr = PetscNewLog(snes,SNESKSPEW,&kctx);CHKERRQ(ierr); 872 snes->kspconvctx = (void*)kctx; 873 kctx->version = 2; 874 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 875 this was too large for some test cases */ 876 kctx->rtol_last = 0; 877 kctx->rtol_max = .9; 878 kctx->gamma = 1.0; 879 kctx->alpha = .5*(1.0 + sqrt(5.0)); 880 kctx->alpha2 = kctx->alpha; 881 kctx->threshold = .1; 882 kctx->lresid_last = 0; 883 kctx->norm_last = 0; 884 885 *outsnes = snes; 886 ierr = PetscPublishAll(snes);CHKERRQ(ierr); 887 PetscFunctionReturn(0); 888 } 889 890 #undef __FUNCT__ 891 #define __FUNCT__ "SNESSetFunction" 892 /*@C 893 SNESSetFunction - Sets the function evaluation routine and function 894 vector for use by the SNES routines in solving systems of nonlinear 895 equations. 896 897 Collective on SNES 898 899 Input Parameters: 900 + snes - the SNES context 901 . r - vector to store function value 902 . func - function evaluation routine 903 - ctx - [optional] user-defined context for private data for the 904 function evaluation routine (may be PETSC_NULL) 905 906 Calling sequence of func: 907 $ func (SNES snes,Vec x,Vec f,void *ctx); 908 909 . f - function vector 910 - ctx - optional user-defined function context 911 912 Notes: 913 The Newton-like methods typically solve linear systems of the form 914 $ f'(x) x = -f(x), 915 where f'(x) denotes the Jacobian matrix and f(x) is the function. 916 917 Level: beginner 918 919 .keywords: SNES, nonlinear, set, function 920 921 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 922 @*/ 923 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 924 { 925 PetscErrorCode ierr; 926 PetscFunctionBegin; 927 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 928 PetscValidHeaderSpecific(r,VEC_COOKIE,2); 929 PetscCheckSameComm(snes,1,r,2); 930 ierr = PetscObjectReference((PetscObject)r);CHKERRQ(ierr); 931 if (snes->vec_func) { ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr); } 932 snes->ops->computefunction = func; 933 snes->vec_func = r; 934 snes->funP = ctx; 935 PetscFunctionReturn(0); 936 } 937 938 /* --------------------------------------------------------------- */ 939 #undef __FUNCT__ 940 #define __FUNCT__ "SNESGetRhs" 941 /*@C 942 SNESGetRhs - Gets the vector for solving F(x) = rhs. If rhs is not set 943 it assumes a zero right hand side. 944 945 Collective on SNES 946 947 Input Parameter: 948 . snes - the SNES context 949 950 Output Parameter: 951 . rhs - the right hand side vector or PETSC_NULL if the right hand side vector is null 952 953 Level: intermediate 954 955 .keywords: SNES, nonlinear, get, function, right hand side 956 957 .seealso: SNESGetSolution(), SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 958 @*/ 959 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetRhs(SNES snes,Vec *rhs) 960 { 961 PetscFunctionBegin; 962 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 963 PetscValidPointer(rhs,2); 964 *rhs = snes->vec_rhs; 965 PetscFunctionReturn(0); 966 } 967 968 #undef __FUNCT__ 969 #define __FUNCT__ "SNESComputeFunction" 970 /*@ 971 SNESComputeFunction - Calls the function that has been set with 972 SNESSetFunction(). 973 974 Collective on SNES 975 976 Input Parameters: 977 + snes - the SNES context 978 - x - input vector 979 980 Output Parameter: 981 . y - function vector, as set by SNESSetFunction() 982 983 Notes: 984 SNESComputeFunction() is typically used within nonlinear solvers 985 implementations, so most users would not generally call this routine 986 themselves. 987 988 Level: developer 989 990 .keywords: SNES, nonlinear, compute, function 991 992 .seealso: SNESSetFunction(), SNESGetFunction() 993 @*/ 994 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeFunction(SNES snes,Vec x,Vec y) 995 { 996 PetscErrorCode ierr; 997 998 PetscFunctionBegin; 999 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1000 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1001 PetscValidHeaderSpecific(y,VEC_COOKIE,3); 1002 PetscCheckSameComm(snes,1,x,2); 1003 PetscCheckSameComm(snes,1,y,3); 1004 1005 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1006 if (snes->ops->computefunction) { 1007 PetscStackPush("SNES user function"); 1008 CHKMEMQ; 1009 ierr = (*snes->ops->computefunction)(snes,x,y,snes->funP); 1010 CHKMEMQ; 1011 PetscStackPop; 1012 if (PetscExceptionValue(ierr)) { 1013 PetscErrorCode pierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(pierr); 1014 } 1015 CHKERRQ(ierr); 1016 } else if (snes->vec_rhs) { 1017 ierr = MatMult(snes->jacobian, x, y);CHKERRQ(ierr); 1018 } else { 1019 SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "Must call SNESSetFunction() before SNESComputeFunction(), likely called from SNESSolve()."); 1020 } 1021 if (snes->vec_rhs) { 1022 ierr = VecAXPY(y,-1.0,snes->vec_rhs);CHKERRQ(ierr); 1023 } 1024 snes->nfuncs++; 1025 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 1026 PetscFunctionReturn(0); 1027 } 1028 1029 #undef __FUNCT__ 1030 #define __FUNCT__ "SNESComputeJacobian" 1031 /*@ 1032 SNESComputeJacobian - Computes the Jacobian matrix that has been 1033 set with SNESSetJacobian(). 1034 1035 Collective on SNES and Mat 1036 1037 Input Parameters: 1038 + snes - the SNES context 1039 - x - input vector 1040 1041 Output Parameters: 1042 + A - Jacobian matrix 1043 . B - optional preconditioning matrix 1044 - flag - flag indicating matrix structure (one of, SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER) 1045 1046 Options Database Keys: 1047 + -snes_lag_preconditioner <lag> 1048 - -snes_lag_jacobian <lag> 1049 1050 Notes: 1051 Most users should not need to explicitly call this routine, as it 1052 is used internally within the nonlinear solvers. 1053 1054 See KSPSetOperators() for important information about setting the 1055 flag parameter. 1056 1057 Level: developer 1058 1059 .keywords: SNES, compute, Jacobian, matrix 1060 1061 .seealso: SNESSetJacobian(), KSPSetOperators(), MatStructure, SNESSetLagPreconditioner(), SNESSetLagJacobian() 1062 @*/ 1063 PetscErrorCode PETSCSNES_DLLEXPORT SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 1064 { 1065 PetscErrorCode ierr; 1066 PetscTruth flag; 1067 1068 PetscFunctionBegin; 1069 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1070 PetscValidHeaderSpecific(X,VEC_COOKIE,2); 1071 PetscValidPointer(flg,5); 1072 PetscCheckSameComm(snes,1,X,2); 1073 if (!snes->ops->computejacobian) PetscFunctionReturn(0); 1074 1075 /* make sure that MatAssemblyBegin/End() is called on A matrix if it is matrix free */ 1076 1077 if (snes->lagjacobian == -2) { 1078 snes->lagjacobian = -1; 1079 ierr = PetscInfo(snes,"Recomputing Jacobian/preconditioner because lag is -2 (means compute Jacobian, but then never again) \n");CHKERRQ(ierr); 1080 } else if (snes->lagjacobian == -1) { 1081 *flg = SAME_PRECONDITIONER; 1082 ierr = PetscInfo(snes,"Reusing Jacobian/preconditioner because lag is -1\n");CHKERRQ(ierr); 1083 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1084 if (flag) { 1085 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1086 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1087 } 1088 PetscFunctionReturn(0); 1089 } else if (snes->lagjacobian > 1 && snes->iter % snes->lagjacobian) { 1090 *flg = SAME_PRECONDITIONER; 1091 ierr = PetscInfo2(snes,"Reusing Jacobian/preconditioner because lag is %D and SNES iteration is %D\n",snes->lagjacobian,snes->iter);CHKERRQ(ierr); 1092 ierr = PetscTypeCompare((PetscObject)*A,MATMFFD,&flag);CHKERRQ(ierr); 1093 if (flag) { 1094 ierr = MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1095 ierr = MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1096 } 1097 PetscFunctionReturn(0); 1098 } 1099 1100 *flg = DIFFERENT_NONZERO_PATTERN; 1101 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1102 PetscStackPush("SNES user Jacobian function"); 1103 CHKMEMQ; 1104 ierr = (*snes->ops->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 1105 CHKMEMQ; 1106 PetscStackPop; 1107 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 1108 1109 if (snes->lagpreconditioner == -2) { 1110 ierr = PetscInfo(snes,"Rebuilding preconditioner exactly once since lag is -2\n");CHKERRQ(ierr); 1111 snes->lagpreconditioner = -1; 1112 } else if (snes->lagpreconditioner == -1) { 1113 *flg = SAME_PRECONDITIONER; 1114 ierr = PetscInfo(snes,"Reusing preconditioner because lag is -1\n");CHKERRQ(ierr); 1115 } else if (snes->lagpreconditioner > 1 && snes->iter % snes->lagpreconditioner) { 1116 *flg = SAME_PRECONDITIONER; 1117 ierr = PetscInfo2(snes,"Reusing preconditioner because lag is %D and SNES iteration is %D\n",snes->lagpreconditioner,snes->iter);CHKERRQ(ierr); 1118 } 1119 1120 /* make sure user returned a correct Jacobian and preconditioner */ 1121 /* PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 1122 PetscValidHeaderSpecific(*B,MAT_COOKIE,4); */ 1123 PetscFunctionReturn(0); 1124 } 1125 1126 #undef __FUNCT__ 1127 #define __FUNCT__ "SNESSetJacobian" 1128 /*@C 1129 SNESSetJacobian - Sets the function to compute Jacobian as well as the 1130 location to store the matrix. 1131 1132 Collective on SNES and Mat 1133 1134 Input Parameters: 1135 + snes - the SNES context 1136 . A - Jacobian matrix 1137 . B - preconditioner matrix (usually same as the Jacobian) 1138 . func - Jacobian evaluation routine 1139 - ctx - [optional] user-defined context for private data for the 1140 Jacobian evaluation routine (may be PETSC_NULL) 1141 1142 Calling sequence of func: 1143 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 1144 1145 + x - input vector 1146 . A - Jacobian matrix 1147 . B - preconditioner matrix, usually the same as A 1148 . flag - flag indicating information about the preconditioner matrix 1149 structure (same as flag in KSPSetOperators()), one of SAME_NONZERO_PATTERN,DIFFERENT_NONZERO_PATTERN,SAME_PRECONDITIONER 1150 - ctx - [optional] user-defined Jacobian context 1151 1152 Notes: 1153 See KSPSetOperators() for important information about setting the flag 1154 output parameter in the routine func(). Be sure to read this information! 1155 1156 The routine func() takes Mat * as the matrix arguments rather than Mat. 1157 This allows the Jacobian evaluation routine to replace A and/or B with a 1158 completely new new matrix structure (not just different matrix elements) 1159 when appropriate, for instance, if the nonzero structure is changing 1160 throughout the global iterations. 1161 1162 If the A matrix and B matrix are different you must call MatAssemblyBegin/End() on 1163 each matrix. 1164 1165 Level: beginner 1166 1167 .keywords: SNES, nonlinear, set, Jacobian, matrix 1168 1169 .seealso: KSPSetOperators(), SNESSetFunction(), MatMFFDComputeJacobian(), SNESDefaultComputeJacobianColor(), MatStructure 1170 @*/ 1171 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 1172 { 1173 PetscErrorCode ierr; 1174 1175 PetscFunctionBegin; 1176 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1177 if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 1178 if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 1179 if (A) PetscCheckSameComm(snes,1,A,2); 1180 if (B) PetscCheckSameComm(snes,1,B,2); 1181 if (func) snes->ops->computejacobian = func; 1182 if (ctx) snes->jacP = ctx; 1183 if (A) { 1184 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 1185 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1186 snes->jacobian = A; 1187 } 1188 if (B) { 1189 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 1190 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1191 snes->jacobian_pre = B; 1192 } 1193 PetscFunctionReturn(0); 1194 } 1195 1196 #undef __FUNCT__ 1197 #define __FUNCT__ "SNESGetJacobian" 1198 /*@C 1199 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 1200 provided context for evaluating the Jacobian. 1201 1202 Not Collective, but Mat object will be parallel if SNES object is 1203 1204 Input Parameter: 1205 . snes - the nonlinear solver context 1206 1207 Output Parameters: 1208 + A - location to stash Jacobian matrix (or PETSC_NULL) 1209 . B - location to stash preconditioner matrix (or PETSC_NULL) 1210 . func - location to put Jacobian function (or PETSC_NULL) 1211 - ctx - location to stash Jacobian ctx (or PETSC_NULL) 1212 1213 Level: advanced 1214 1215 .seealso: SNESSetJacobian(), SNESComputeJacobian() 1216 @*/ 1217 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetJacobian(SNES snes,Mat *A,Mat *B,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void **ctx) 1218 { 1219 PetscFunctionBegin; 1220 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1221 if (A) *A = snes->jacobian; 1222 if (B) *B = snes->jacobian_pre; 1223 if (func) *func = snes->ops->computejacobian; 1224 if (ctx) *ctx = snes->jacP; 1225 PetscFunctionReturn(0); 1226 } 1227 1228 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 1229 EXTERN PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 1230 1231 #undef __FUNCT__ 1232 #define __FUNCT__ "SNESSetUp" 1233 /*@ 1234 SNESSetUp - Sets up the internal data structures for the later use 1235 of a nonlinear solver. 1236 1237 Collective on SNES 1238 1239 Input Parameters: 1240 . snes - the SNES context 1241 1242 Notes: 1243 For basic use of the SNES solvers the user need not explicitly call 1244 SNESSetUp(), since these actions will automatically occur during 1245 the call to SNESSolve(). However, if one wishes to control this 1246 phase separately, SNESSetUp() should be called after SNESCreate() 1247 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 1248 1249 Level: advanced 1250 1251 .keywords: SNES, nonlinear, setup 1252 1253 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 1254 @*/ 1255 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUp(SNES snes) 1256 { 1257 PetscErrorCode ierr; 1258 PetscTruth flg; 1259 1260 PetscFunctionBegin; 1261 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1262 if (snes->setupcalled) PetscFunctionReturn(0); 1263 1264 if (!((PetscObject)snes)->type_name) { 1265 ierr = SNESSetType(snes,SNESLS);CHKERRQ(ierr); 1266 } 1267 1268 ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 1269 /* 1270 This version replaces the user provided Jacobian matrix with a 1271 matrix-free version but still employs the user-provided preconditioner matrix 1272 */ 1273 if (flg) { 1274 Mat J; 1275 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 1276 ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr); 1277 ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr); 1278 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1279 ierr = MatDestroy(J);CHKERRQ(ierr); 1280 } 1281 1282 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) && !defined(PETSC_USE_MAT_SINGLE) && !defined(PETSC_USE_LONG_DOUBLE) && !defined(PETSC_USE_INT) 1283 ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 1284 if (flg) { 1285 Mat J; 1286 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1287 ierr = PetscInfo(snes,"Setting default matrix-free operator routines (version 2)\n");CHKERRQ(ierr); 1288 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 1289 ierr = MatDestroy(J);CHKERRQ(ierr); 1290 } 1291 #endif 1292 1293 ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 1294 /* 1295 This version replaces both the user-provided Jacobian and the user- 1296 provided preconditioner matrix with the default matrix free version. 1297 */ 1298 if (flg) { 1299 Mat J; 1300 KSP ksp; 1301 PC pc; 1302 /* create and set matrix-free operator */ 1303 ierr = MatCreateSNESMF(snes,&J);CHKERRQ(ierr); 1304 ierr = MatMFFDSetFromOptions(J);CHKERRQ(ierr); 1305 ierr = PetscInfo(snes,"Setting default matrix-free operator routines\n");CHKERRQ(ierr); 1306 ierr = SNESSetJacobian(snes,J,J,MatMFFDComputeJacobian,snes->funP);CHKERRQ(ierr); 1307 ierr = MatDestroy(J);CHKERRQ(ierr); 1308 /* force no preconditioner */ 1309 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1310 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1311 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1312 if (!flg) { 1313 ierr = PetscInfo(snes,"Setting default matrix-free preconditioner routines;\nThat is no preconditioner is being used\n");CHKERRQ(ierr); 1314 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1315 } 1316 } 1317 1318 if (!snes->vec_func && !snes->vec_rhs) { 1319 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1320 } 1321 if (!snes->ops->computefunction && !snes->vec_rhs) { 1322 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1323 } 1324 if (!snes->jacobian) { 1325 SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1326 } 1327 if (snes->vec_func == snes->vec_sol) { 1328 SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1329 } 1330 1331 if (snes->ops->setup) { 1332 ierr = (*snes->ops->setup)(snes);CHKERRQ(ierr); 1333 } 1334 snes->setupcalled = PETSC_TRUE; 1335 PetscFunctionReturn(0); 1336 } 1337 1338 #undef __FUNCT__ 1339 #define __FUNCT__ "SNESDestroy" 1340 /*@ 1341 SNESDestroy - Destroys the nonlinear solver context that was created 1342 with SNESCreate(). 1343 1344 Collective on SNES 1345 1346 Input Parameter: 1347 . snes - the SNES context 1348 1349 Level: beginner 1350 1351 .keywords: SNES, nonlinear, destroy 1352 1353 .seealso: SNESCreate(), SNESSolve() 1354 @*/ 1355 PetscErrorCode PETSCSNES_DLLEXPORT SNESDestroy(SNES snes) 1356 { 1357 PetscErrorCode ierr; 1358 1359 PetscFunctionBegin; 1360 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1361 if (--((PetscObject)snes)->refct > 0) PetscFunctionReturn(0); 1362 1363 /* if memory was published with AMS then destroy it */ 1364 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1365 if (snes->ops->destroy) {ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr);} 1366 1367 if (snes->vec_rhs) {ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr);} 1368 if (snes->vec_sol) {ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr);} 1369 if (snes->vec_func) {ierr = VecDestroy(snes->vec_func);CHKERRQ(ierr);} 1370 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1371 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1372 if (snes->ksp) {ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr);} 1373 ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr); 1374 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1375 ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 1376 if (snes->ops->convergeddestroy) {ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr);} 1377 ierr = PetscHeaderDestroy(snes);CHKERRQ(ierr); 1378 PetscFunctionReturn(0); 1379 } 1380 1381 /* ----------- Routines to set solver parameters ---------- */ 1382 1383 #undef __FUNCT__ 1384 #define __FUNCT__ "SNESSetLagPreconditioner" 1385 /*@ 1386 SNESSetLagPreconditioner - Determines when the preconditioner is rebuilt in the nonlinear solve. 1387 1388 Collective on SNES 1389 1390 Input Parameters: 1391 + snes - the SNES context 1392 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1393 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 1394 1395 Options Database Keys: 1396 . -snes_lag_preconditioner <lag> 1397 1398 Notes: 1399 The default is 1 1400 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1401 If -1 is used before the very first nonlinear solve the preconditioner is still built because there is no previous preconditioner to use 1402 1403 Level: intermediate 1404 1405 .keywords: SNES, nonlinear, set, convergence, tolerances 1406 1407 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagJacobian(), SNESGetLagJacobian() 1408 1409 @*/ 1410 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagPreconditioner(SNES snes,PetscInt lag) 1411 { 1412 PetscFunctionBegin; 1413 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1414 if (lag < -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 1415 if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1416 snes->lagpreconditioner = lag; 1417 PetscFunctionReturn(0); 1418 } 1419 1420 #undef __FUNCT__ 1421 #define __FUNCT__ "SNESGetLagPreconditioner" 1422 /*@ 1423 SNESGetLagPreconditioner - Indicates how often the preconditioner is rebuilt 1424 1425 Collective on SNES 1426 1427 Input Parameter: 1428 . snes - the SNES context 1429 1430 Output Parameter: 1431 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1432 the Jacobian is built etc. -2 indicates rebuild preconditioner at next chance but then never rebuild after that 1433 1434 Options Database Keys: 1435 . -snes_lag_preconditioner <lag> 1436 1437 Notes: 1438 The default is 1 1439 The preconditioner is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1440 1441 Level: intermediate 1442 1443 .keywords: SNES, nonlinear, set, convergence, tolerances 1444 1445 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagPreconditioner() 1446 1447 @*/ 1448 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagPreconditioner(SNES snes,PetscInt *lag) 1449 { 1450 PetscFunctionBegin; 1451 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1452 *lag = snes->lagpreconditioner; 1453 PetscFunctionReturn(0); 1454 } 1455 1456 #undef __FUNCT__ 1457 #define __FUNCT__ "SNESSetLagJacobian" 1458 /*@ 1459 SNESSetLagJacobian - Determines when the Jacobian is rebuilt in the nonlinear solve. See SNESSetLagPreconditioner() for determining how 1460 often the preconditioner is rebuilt. 1461 1462 Collective on SNES 1463 1464 Input Parameters: 1465 + snes - the SNES context 1466 - lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1467 the Jacobian is built etc. -2 means rebuild at next chance but then never again 1468 1469 Options Database Keys: 1470 . -snes_lag_jacobian <lag> 1471 1472 Notes: 1473 The default is 1 1474 The Jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1475 If -1 is used before the very first nonlinear solve the CODE WILL FAIL! because no Jacobian is used, use -2 to indicate you want it recomputed 1476 at the next Newton step but never again (unless it is reset to another value) 1477 1478 Level: intermediate 1479 1480 .keywords: SNES, nonlinear, set, convergence, tolerances 1481 1482 .seealso: SNESSetTrustRegionTolerance(), SNESGetLagPreconditioner(), SNESSetLagPreconditioner(), SNESGetLagJacobian() 1483 1484 @*/ 1485 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetLagJacobian(SNES snes,PetscInt lag) 1486 { 1487 PetscFunctionBegin; 1488 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1489 if (lag < -2) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag must be -2, -1, 1 or greater"); 1490 if (!lag) SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Lag cannot be 0"); 1491 snes->lagjacobian = lag; 1492 PetscFunctionReturn(0); 1493 } 1494 1495 #undef __FUNCT__ 1496 #define __FUNCT__ "SNESGetLagJacobian" 1497 /*@ 1498 SNESGetLagJacobian - Indicates how often the Jacobian is rebuilt. See SNESGetLagPreconditioner() to determine when the preconditioner is rebuilt 1499 1500 Collective on SNES 1501 1502 Input Parameter: 1503 . snes - the SNES context 1504 1505 Output Parameter: 1506 . lag - -1 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is computed within a single nonlinear solve, 2 means every second time 1507 the Jacobian is built etc. 1508 1509 Options Database Keys: 1510 . -snes_lag_jacobian <lag> 1511 1512 Notes: 1513 The default is 1 1514 The jacobian is ALWAYS built in the first iteration of a nonlinear solve unless lag is -1 1515 1516 Level: intermediate 1517 1518 .keywords: SNES, nonlinear, set, convergence, tolerances 1519 1520 .seealso: SNESSetTrustRegionTolerance(), SNESSetLagJacobian(), SNESSetLagPreconditioner(), SNESGetLagPreconditioner() 1521 1522 @*/ 1523 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetLagJacobian(SNES snes,PetscInt *lag) 1524 { 1525 PetscFunctionBegin; 1526 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1527 *lag = snes->lagjacobian; 1528 PetscFunctionReturn(0); 1529 } 1530 1531 #undef __FUNCT__ 1532 #define __FUNCT__ "SNESSetTolerances" 1533 /*@ 1534 SNESSetTolerances - Sets various parameters used in convergence tests. 1535 1536 Collective on SNES 1537 1538 Input Parameters: 1539 + snes - the SNES context 1540 . abstol - absolute convergence tolerance 1541 . rtol - relative convergence tolerance 1542 . stol - convergence tolerance in terms of the norm 1543 of the change in the solution between steps 1544 . maxit - maximum number of iterations 1545 - maxf - maximum number of function evaluations 1546 1547 Options Database Keys: 1548 + -snes_atol <abstol> - Sets abstol 1549 . -snes_rtol <rtol> - Sets rtol 1550 . -snes_stol <stol> - Sets stol 1551 . -snes_max_it <maxit> - Sets maxit 1552 - -snes_max_funcs <maxf> - Sets maxf 1553 1554 Notes: 1555 The default maximum number of iterations is 50. 1556 The default maximum number of function evaluations is 1000. 1557 1558 Level: intermediate 1559 1560 .keywords: SNES, nonlinear, set, convergence, tolerances 1561 1562 .seealso: SNESSetTrustRegionTolerance() 1563 @*/ 1564 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1565 { 1566 PetscFunctionBegin; 1567 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1568 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1569 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1570 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1571 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1572 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1573 PetscFunctionReturn(0); 1574 } 1575 1576 #undef __FUNCT__ 1577 #define __FUNCT__ "SNESGetTolerances" 1578 /*@ 1579 SNESGetTolerances - Gets various parameters used in convergence tests. 1580 1581 Not Collective 1582 1583 Input Parameters: 1584 + snes - the SNES context 1585 . atol - absolute convergence tolerance 1586 . rtol - relative convergence tolerance 1587 . stol - convergence tolerance in terms of the norm 1588 of the change in the solution between steps 1589 . maxit - maximum number of iterations 1590 - maxf - maximum number of function evaluations 1591 1592 Notes: 1593 The user can specify PETSC_NULL for any parameter that is not needed. 1594 1595 Level: intermediate 1596 1597 .keywords: SNES, nonlinear, get, convergence, tolerances 1598 1599 .seealso: SNESSetTolerances() 1600 @*/ 1601 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1602 { 1603 PetscFunctionBegin; 1604 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1605 if (atol) *atol = snes->abstol; 1606 if (rtol) *rtol = snes->rtol; 1607 if (stol) *stol = snes->xtol; 1608 if (maxit) *maxit = snes->max_its; 1609 if (maxf) *maxf = snes->max_funcs; 1610 PetscFunctionReturn(0); 1611 } 1612 1613 #undef __FUNCT__ 1614 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1615 /*@ 1616 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1617 1618 Collective on SNES 1619 1620 Input Parameters: 1621 + snes - the SNES context 1622 - tol - tolerance 1623 1624 Options Database Key: 1625 . -snes_trtol <tol> - Sets tol 1626 1627 Level: intermediate 1628 1629 .keywords: SNES, nonlinear, set, trust region, tolerance 1630 1631 .seealso: SNESSetTolerances() 1632 @*/ 1633 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1634 { 1635 PetscFunctionBegin; 1636 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1637 snes->deltatol = tol; 1638 PetscFunctionReturn(0); 1639 } 1640 1641 /* 1642 Duplicate the lg monitors for SNES from KSP; for some reason with 1643 dynamic libraries things don't work under Sun4 if we just use 1644 macros instead of functions 1645 */ 1646 #undef __FUNCT__ 1647 #define __FUNCT__ "SNESMonitorLG" 1648 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLG(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1649 { 1650 PetscErrorCode ierr; 1651 1652 PetscFunctionBegin; 1653 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1654 ierr = KSPMonitorLG((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1655 PetscFunctionReturn(0); 1656 } 1657 1658 #undef __FUNCT__ 1659 #define __FUNCT__ "SNESMonitorLGCreate" 1660 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1661 { 1662 PetscErrorCode ierr; 1663 1664 PetscFunctionBegin; 1665 ierr = KSPMonitorLGCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1666 PetscFunctionReturn(0); 1667 } 1668 1669 #undef __FUNCT__ 1670 #define __FUNCT__ "SNESMonitorLGDestroy" 1671 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorLGDestroy(PetscDrawLG draw) 1672 { 1673 PetscErrorCode ierr; 1674 1675 PetscFunctionBegin; 1676 ierr = KSPMonitorLGDestroy(draw);CHKERRQ(ierr); 1677 PetscFunctionReturn(0); 1678 } 1679 1680 /* ------------ Routines to set performance monitoring options ----------- */ 1681 1682 #undef __FUNCT__ 1683 #define __FUNCT__ "SNESMonitorSet" 1684 /*@C 1685 SNESMonitorSet - Sets an ADDITIONAL function that is to be used at every 1686 iteration of the nonlinear solver to display the iteration's 1687 progress. 1688 1689 Collective on SNES 1690 1691 Input Parameters: 1692 + snes - the SNES context 1693 . func - monitoring routine 1694 . mctx - [optional] user-defined context for private data for the 1695 monitor routine (use PETSC_NULL if no context is desired) 1696 - monitordestroy - [optional] routine that frees monitor context 1697 (may be PETSC_NULL) 1698 1699 Calling sequence of func: 1700 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1701 1702 + snes - the SNES context 1703 . its - iteration number 1704 . norm - 2-norm function value (may be estimated) 1705 - mctx - [optional] monitoring context 1706 1707 Options Database Keys: 1708 + -snes_monitor - sets SNESMonitorDefault() 1709 . -snes_monitor_draw - sets line graph monitor, 1710 uses SNESMonitorLGCreate() 1711 _ -snes_monitor_cancel - cancels all monitors that have 1712 been hardwired into a code by 1713 calls to SNESMonitorSet(), but 1714 does not cancel those set via 1715 the options database. 1716 1717 Notes: 1718 Several different monitoring routines may be set by calling 1719 SNESMonitorSet() multiple times; all will be called in the 1720 order in which they were set. 1721 1722 Fortran notes: Only a single monitor function can be set for each SNES object 1723 1724 Level: intermediate 1725 1726 .keywords: SNES, nonlinear, set, monitor 1727 1728 .seealso: SNESMonitorDefault(), SNESMonitorCancel() 1729 @*/ 1730 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorSet(SNES snes,PetscErrorCode (*monitor)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 1731 { 1732 PetscInt i; 1733 1734 PetscFunctionBegin; 1735 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1736 if (snes->numbermonitors >= MAXSNESMONITORS) { 1737 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1738 } 1739 for (i=0; i<snes->numbermonitors;i++) { 1740 if (monitor == snes->monitor[i] && monitordestroy == snes->monitordestroy[i] && mctx == snes->monitorcontext[i]) PetscFunctionReturn(0); 1741 1742 /* check if both default monitors that share common ASCII viewer */ 1743 if (monitor == snes->monitor[i] && monitor == SNESMonitorDefault) { 1744 if (mctx && snes->monitorcontext[i]) { 1745 PetscErrorCode ierr; 1746 PetscViewerASCIIMonitor viewer1 = (PetscViewerASCIIMonitor) mctx; 1747 PetscViewerASCIIMonitor viewer2 = (PetscViewerASCIIMonitor) snes->monitorcontext[i]; 1748 if (viewer1->viewer == viewer2->viewer) { 1749 ierr = (*monitordestroy)(mctx);CHKERRQ(ierr); 1750 PetscFunctionReturn(0); 1751 } 1752 } 1753 } 1754 } 1755 snes->monitor[snes->numbermonitors] = monitor; 1756 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1757 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1758 PetscFunctionReturn(0); 1759 } 1760 1761 #undef __FUNCT__ 1762 #define __FUNCT__ "SNESMonitorCancel" 1763 /*@C 1764 SNESMonitorCancel - Clears all the monitor functions for a SNES object. 1765 1766 Collective on SNES 1767 1768 Input Parameters: 1769 . snes - the SNES context 1770 1771 Options Database Key: 1772 . -snes_monitor_cancel - cancels all monitors that have been hardwired 1773 into a code by calls to SNESMonitorSet(), but does not cancel those 1774 set via the options database 1775 1776 Notes: 1777 There is no way to clear one specific monitor from a SNES object. 1778 1779 Level: intermediate 1780 1781 .keywords: SNES, nonlinear, set, monitor 1782 1783 .seealso: SNESMonitorDefault(), SNESMonitorSet() 1784 @*/ 1785 PetscErrorCode PETSCSNES_DLLEXPORT SNESMonitorCancel(SNES snes) 1786 { 1787 PetscErrorCode ierr; 1788 PetscInt i; 1789 1790 PetscFunctionBegin; 1791 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1792 for (i=0; i<snes->numbermonitors; i++) { 1793 if (snes->monitordestroy[i]) { 1794 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1795 } 1796 } 1797 snes->numbermonitors = 0; 1798 PetscFunctionReturn(0); 1799 } 1800 1801 #undef __FUNCT__ 1802 #define __FUNCT__ "SNESSetConvergenceTest" 1803 /*@C 1804 SNESSetConvergenceTest - Sets the function that is to be used 1805 to test for convergence of the nonlinear iterative solution. 1806 1807 Collective on SNES 1808 1809 Input Parameters: 1810 + snes - the SNES context 1811 . func - routine to test for convergence 1812 . cctx - [optional] context for private data for the convergence routine (may be PETSC_NULL) 1813 - destroy - [optional] destructor for the context (may be PETSC_NULL; PETSC_NULL_FUNCTION in Fortran) 1814 1815 Calling sequence of func: 1816 $ PetscErrorCode func (SNES snes,PetscInt it,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1817 1818 + snes - the SNES context 1819 . it - current iteration (0 is the first and is before any Newton step) 1820 . cctx - [optional] convergence context 1821 . reason - reason for convergence/divergence 1822 . xnorm - 2-norm of current iterate 1823 . gnorm - 2-norm of current step 1824 - f - 2-norm of function 1825 1826 Level: advanced 1827 1828 .keywords: SNES, nonlinear, set, convergence, test 1829 1830 .seealso: SNESDefaultConverged(), SNESSkipConverged() 1831 @*/ 1832 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx,PetscErrorCode (*destroy)(void*)) 1833 { 1834 PetscErrorCode ierr; 1835 1836 PetscFunctionBegin; 1837 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1838 if (!func) func = SNESSkipConverged; 1839 if (snes->ops->convergeddestroy) { 1840 ierr = (*snes->ops->convergeddestroy)(snes->cnvP);CHKERRQ(ierr); 1841 } 1842 snes->ops->converged = func; 1843 snes->ops->convergeddestroy = destroy; 1844 snes->cnvP = cctx; 1845 PetscFunctionReturn(0); 1846 } 1847 1848 #undef __FUNCT__ 1849 #define __FUNCT__ "SNESGetConvergedReason" 1850 /*@ 1851 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1852 1853 Not Collective 1854 1855 Input Parameter: 1856 . snes - the SNES context 1857 1858 Output Parameter: 1859 . reason - negative value indicates diverged, positive value converged, see SNESConvergedReason or the 1860 manual pages for the individual convergence tests for complete lists 1861 1862 Level: intermediate 1863 1864 Notes: Can only be called after the call the SNESSolve() is complete. 1865 1866 .keywords: SNES, nonlinear, set, convergence, test 1867 1868 .seealso: SNESSetConvergenceTest(), SNESConvergedReason 1869 @*/ 1870 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1871 { 1872 PetscFunctionBegin; 1873 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1874 PetscValidPointer(reason,2); 1875 *reason = snes->reason; 1876 PetscFunctionReturn(0); 1877 } 1878 1879 #undef __FUNCT__ 1880 #define __FUNCT__ "SNESSetConvergenceHistory" 1881 /*@ 1882 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1883 1884 Collective on SNES 1885 1886 Input Parameters: 1887 + snes - iterative context obtained from SNESCreate() 1888 . a - array to hold history 1889 . its - integer array holds the number of linear iterations for each solve. 1890 . na - size of a and its 1891 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1892 else it continues storing new values for new nonlinear solves after the old ones 1893 1894 Notes: 1895 If set, this array will contain the function norms computed 1896 at each step. 1897 1898 This routine is useful, e.g., when running a code for purposes 1899 of accurate performance monitoring, when no I/O should be done 1900 during the section of code that is being timed. 1901 1902 Level: intermediate 1903 1904 .keywords: SNES, set, convergence, history 1905 1906 .seealso: SNESGetConvergenceHistory() 1907 1908 @*/ 1909 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt its[],PetscInt na,PetscTruth reset) 1910 { 1911 PetscFunctionBegin; 1912 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1913 if (na) PetscValidScalarPointer(a,2); 1914 if (its) PetscValidIntPointer(its,3); 1915 snes->conv_hist = a; 1916 snes->conv_hist_its = its; 1917 snes->conv_hist_max = na; 1918 snes->conv_hist_len = 0; 1919 snes->conv_hist_reset = reset; 1920 PetscFunctionReturn(0); 1921 } 1922 1923 #undef __FUNCT__ 1924 #define __FUNCT__ "SNESGetConvergenceHistory" 1925 /*@C 1926 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1927 1928 Collective on SNES 1929 1930 Input Parameter: 1931 . snes - iterative context obtained from SNESCreate() 1932 1933 Output Parameters: 1934 . a - array to hold history 1935 . its - integer array holds the number of linear iterations (or 1936 negative if not converged) for each solve. 1937 - na - size of a and its 1938 1939 Notes: 1940 The calling sequence for this routine in Fortran is 1941 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1942 1943 This routine is useful, e.g., when running a code for purposes 1944 of accurate performance monitoring, when no I/O should be done 1945 during the section of code that is being timed. 1946 1947 Level: intermediate 1948 1949 .keywords: SNES, get, convergence, history 1950 1951 .seealso: SNESSetConvergencHistory() 1952 1953 @*/ 1954 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1955 { 1956 PetscFunctionBegin; 1957 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1958 if (a) *a = snes->conv_hist; 1959 if (its) *its = snes->conv_hist_its; 1960 if (na) *na = snes->conv_hist_len; 1961 PetscFunctionReturn(0); 1962 } 1963 1964 #undef __FUNCT__ 1965 #define __FUNCT__ "SNESSetUpdate" 1966 /*@C 1967 SNESSetUpdate - Sets the general-purpose update function called 1968 at the beginning o every iteration of the nonlinear solve. Specifically 1969 it is called just before the Jacobian is "evaluated". 1970 1971 Collective on SNES 1972 1973 Input Parameters: 1974 . snes - The nonlinear solver context 1975 . func - The function 1976 1977 Calling sequence of func: 1978 . func (SNES snes, PetscInt step); 1979 1980 . step - The current step of the iteration 1981 1982 Level: intermediate 1983 1984 .keywords: SNES, update 1985 1986 .seealso SNESDefaultUpdate(), SNESSetJacobian(), SNESSolve() 1987 @*/ 1988 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 1989 { 1990 PetscFunctionBegin; 1991 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1992 snes->ops->update = func; 1993 PetscFunctionReturn(0); 1994 } 1995 1996 #undef __FUNCT__ 1997 #define __FUNCT__ "SNESDefaultUpdate" 1998 /*@ 1999 SNESDefaultUpdate - The default update function which does nothing. 2000 2001 Not collective 2002 2003 Input Parameters: 2004 . snes - The nonlinear solver context 2005 . step - The current step of the iteration 2006 2007 Level: intermediate 2008 2009 .keywords: SNES, update 2010 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultShortolutionBC() 2011 @*/ 2012 PetscErrorCode PETSCSNES_DLLEXPORT SNESDefaultUpdate(SNES snes, PetscInt step) 2013 { 2014 PetscFunctionBegin; 2015 PetscFunctionReturn(0); 2016 } 2017 2018 #undef __FUNCT__ 2019 #define __FUNCT__ "SNESScaleStep_Private" 2020 /* 2021 SNESScaleStep_Private - Scales a step so that its length is less than the 2022 positive parameter delta. 2023 2024 Input Parameters: 2025 + snes - the SNES context 2026 . y - approximate solution of linear system 2027 . fnorm - 2-norm of current function 2028 - delta - trust region size 2029 2030 Output Parameters: 2031 + gpnorm - predicted function norm at the new point, assuming local 2032 linearization. The value is zero if the step lies within the trust 2033 region, and exceeds zero otherwise. 2034 - ynorm - 2-norm of the step 2035 2036 Note: 2037 For non-trust region methods such as SNESLS, the parameter delta 2038 is set to be the maximum allowable step size. 2039 2040 .keywords: SNES, nonlinear, scale, step 2041 */ 2042 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 2043 { 2044 PetscReal nrm; 2045 PetscScalar cnorm; 2046 PetscErrorCode ierr; 2047 2048 PetscFunctionBegin; 2049 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2050 PetscValidHeaderSpecific(y,VEC_COOKIE,2); 2051 PetscCheckSameComm(snes,1,y,2); 2052 2053 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 2054 if (nrm > *delta) { 2055 nrm = *delta/nrm; 2056 *gpnorm = (1.0 - nrm)*(*fnorm); 2057 cnorm = nrm; 2058 ierr = VecScale(y,cnorm);CHKERRQ(ierr); 2059 *ynorm = *delta; 2060 } else { 2061 *gpnorm = 0.0; 2062 *ynorm = nrm; 2063 } 2064 PetscFunctionReturn(0); 2065 } 2066 2067 #undef __FUNCT__ 2068 #define __FUNCT__ "SNESSolve" 2069 /*@C 2070 SNESSolve - Solves a nonlinear system F(x) = b. 2071 Call SNESSolve() after calling SNESCreate() and optional routines of the form SNESSetXXX(). 2072 2073 Collective on SNES 2074 2075 Input Parameters: 2076 + snes - the SNES context 2077 . b - the constant part of the equation, or PETSC_NULL to use zero. 2078 - x - the solution vector. 2079 2080 Notes: 2081 The user should initialize the vector,x, with the initial guess 2082 for the nonlinear solve prior to calling SNESSolve. In particular, 2083 to employ an initial guess of zero, the user should explicitly set 2084 this vector to zero by calling VecSet(). 2085 2086 Level: beginner 2087 2088 .keywords: SNES, nonlinear, solve 2089 2090 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian() 2091 @*/ 2092 PetscErrorCode PETSCSNES_DLLEXPORT SNESSolve(SNES snes,Vec b,Vec x) 2093 { 2094 PetscErrorCode ierr; 2095 PetscTruth flg; 2096 char filename[PETSC_MAX_PATH_LEN]; 2097 PetscViewer viewer; 2098 2099 PetscFunctionBegin; 2100 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2101 PetscValidHeaderSpecific(x,VEC_COOKIE,3); 2102 PetscCheckSameComm(snes,1,x,3); 2103 if (b) PetscValidHeaderSpecific(b,VEC_COOKIE,2); 2104 if (b) PetscCheckSameComm(snes,1,b,2); 2105 2106 /* set solution vector */ 2107 ierr = PetscObjectReference((PetscObject)x);CHKERRQ(ierr); 2108 if (snes->vec_sol) { ierr = VecDestroy(snes->vec_sol);CHKERRQ(ierr); } 2109 snes->vec_sol = x; 2110 /* set afine vector if provided */ 2111 if (b) { ierr = PetscObjectReference((PetscObject)b);CHKERRQ(ierr); } 2112 if (snes->vec_rhs) { ierr = VecDestroy(snes->vec_rhs);CHKERRQ(ierr); } 2113 snes->vec_rhs = b; 2114 2115 if (!snes->vec_func && snes->vec_rhs) { 2116 ierr = VecDuplicate(b, &snes->vec_func); CHKERRQ(ierr); 2117 } 2118 2119 ierr = SNESSetUp(snes);CHKERRQ(ierr); 2120 2121 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 2122 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 2123 2124 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2125 ierr = (*snes->ops->solve)(snes);CHKERRQ(ierr); 2126 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 2127 if (snes->domainerror){ 2128 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 2129 snes->domainerror = PETSC_FALSE; 2130 } 2131 2132 if (!snes->reason) { 2133 SETERRQ(PETSC_ERR_PLIB,"Internal error, solver returned without setting converged reason"); 2134 } 2135 2136 ierr = PetscOptionsGetString(((PetscObject)snes)->prefix,"-snes_view",filename,PETSC_MAX_PATH_LEN,&flg);CHKERRQ(ierr); 2137 if (flg && !PetscPreLoadingOn) { 2138 ierr = PetscViewerASCIIOpen(((PetscObject)snes)->comm,filename,&viewer);CHKERRQ(ierr); 2139 ierr = SNESView(snes,viewer);CHKERRQ(ierr); 2140 ierr = PetscViewerDestroy(viewer);CHKERRQ(ierr); 2141 } 2142 2143 ierr = PetscOptionsHasName(((PetscObject)snes)->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 2144 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 2145 if (snes->printreason) { 2146 if (snes->reason > 0) { 2147 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve converged due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2148 } else { 2149 ierr = PetscPrintf(((PetscObject)snes)->comm,"Nonlinear solve did not converge due to %s\n",SNESConvergedReasons[snes->reason]);CHKERRQ(ierr); 2150 } 2151 } 2152 2153 PetscFunctionReturn(0); 2154 } 2155 2156 /* --------- Internal routines for SNES Package --------- */ 2157 2158 #undef __FUNCT__ 2159 #define __FUNCT__ "SNESSetType" 2160 /*@C 2161 SNESSetType - Sets the method for the nonlinear solver. 2162 2163 Collective on SNES 2164 2165 Input Parameters: 2166 + snes - the SNES context 2167 - type - a known method 2168 2169 Options Database Key: 2170 . -snes_type <type> - Sets the method; use -help for a list 2171 of available methods (for instance, ls or tr) 2172 2173 Notes: 2174 See "petsc/include/petscsnes.h" for available methods (for instance) 2175 + SNESLS - Newton's method with line search 2176 (systems of nonlinear equations) 2177 . SNESTR - Newton's method with trust region 2178 (systems of nonlinear equations) 2179 2180 Normally, it is best to use the SNESSetFromOptions() command and then 2181 set the SNES solver type from the options database rather than by using 2182 this routine. Using the options database provides the user with 2183 maximum flexibility in evaluating the many nonlinear solvers. 2184 The SNESSetType() routine is provided for those situations where it 2185 is necessary to set the nonlinear solver independently of the command 2186 line or options database. This might be the case, for example, when 2187 the choice of solver changes during the execution of the program, 2188 and the user's application is taking responsibility for choosing the 2189 appropriate method. 2190 2191 Level: intermediate 2192 2193 .keywords: SNES, set, type 2194 2195 .seealso: SNESType, SNESCreate() 2196 2197 @*/ 2198 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetType(SNES snes,const SNESType type) 2199 { 2200 PetscErrorCode ierr,(*r)(SNES); 2201 PetscTruth match; 2202 2203 PetscFunctionBegin; 2204 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2205 PetscValidCharPointer(type,2); 2206 2207 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 2208 if (match) PetscFunctionReturn(0); 2209 2210 ierr = PetscFListFind(SNESList,((PetscObject)snes)->comm,type,(void (**)(void)) &r);CHKERRQ(ierr); 2211 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 2212 /* Destroy the previous private SNES context */ 2213 if (snes->ops->destroy) { ierr = (*(snes)->ops->destroy)(snes);CHKERRQ(ierr); } 2214 /* Reinitialize function pointers in SNESOps structure */ 2215 snes->ops->setup = 0; 2216 snes->ops->solve = 0; 2217 snes->ops->view = 0; 2218 snes->ops->setfromoptions = 0; 2219 snes->ops->destroy = 0; 2220 /* Call the SNESCreate_XXX routine for this particular Nonlinear solver */ 2221 snes->setupcalled = PETSC_FALSE; 2222 ierr = (*r)(snes);CHKERRQ(ierr); 2223 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 2224 PetscFunctionReturn(0); 2225 } 2226 2227 2228 /* --------------------------------------------------------------------- */ 2229 #undef __FUNCT__ 2230 #define __FUNCT__ "SNESRegisterDestroy" 2231 /*@ 2232 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 2233 registered by SNESRegisterDynamic(). 2234 2235 Not Collective 2236 2237 Level: advanced 2238 2239 .keywords: SNES, nonlinear, register, destroy 2240 2241 .seealso: SNESRegisterAll(), SNESRegisterAll() 2242 @*/ 2243 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegisterDestroy(void) 2244 { 2245 PetscErrorCode ierr; 2246 2247 PetscFunctionBegin; 2248 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 2249 SNESRegisterAllCalled = PETSC_FALSE; 2250 PetscFunctionReturn(0); 2251 } 2252 2253 #undef __FUNCT__ 2254 #define __FUNCT__ "SNESGetType" 2255 /*@C 2256 SNESGetType - Gets the SNES method type and name (as a string). 2257 2258 Not Collective 2259 2260 Input Parameter: 2261 . snes - nonlinear solver context 2262 2263 Output Parameter: 2264 . type - SNES method (a character string) 2265 2266 Level: intermediate 2267 2268 .keywords: SNES, nonlinear, get, type, name 2269 @*/ 2270 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetType(SNES snes,const SNESType *type) 2271 { 2272 PetscFunctionBegin; 2273 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2274 PetscValidPointer(type,2); 2275 *type = ((PetscObject)snes)->type_name; 2276 PetscFunctionReturn(0); 2277 } 2278 2279 #undef __FUNCT__ 2280 #define __FUNCT__ "SNESGetSolution" 2281 /*@ 2282 SNESGetSolution - Returns the vector where the approximate solution is 2283 stored. 2284 2285 Not Collective, but Vec is parallel if SNES is parallel 2286 2287 Input Parameter: 2288 . snes - the SNES context 2289 2290 Output Parameter: 2291 . x - the solution 2292 2293 Level: intermediate 2294 2295 .keywords: SNES, nonlinear, get, solution 2296 2297 .seealso: SNESGetSolutionUpdate(), SNESGetFunction() 2298 @*/ 2299 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolution(SNES snes,Vec *x) 2300 { 2301 PetscFunctionBegin; 2302 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2303 PetscValidPointer(x,2); 2304 *x = snes->vec_sol; 2305 PetscFunctionReturn(0); 2306 } 2307 2308 #undef __FUNCT__ 2309 #define __FUNCT__ "SNESGetSolutionUpdate" 2310 /*@ 2311 SNESGetSolutionUpdate - Returns the vector where the solution update is 2312 stored. 2313 2314 Not Collective, but Vec is parallel if SNES is parallel 2315 2316 Input Parameter: 2317 . snes - the SNES context 2318 2319 Output Parameter: 2320 . x - the solution update 2321 2322 Level: advanced 2323 2324 .keywords: SNES, nonlinear, get, solution, update 2325 2326 .seealso: SNESGetSolution(), SNESGetFunction() 2327 @*/ 2328 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetSolutionUpdate(SNES snes,Vec *x) 2329 { 2330 PetscFunctionBegin; 2331 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2332 PetscValidPointer(x,2); 2333 *x = snes->vec_sol_update; 2334 PetscFunctionReturn(0); 2335 } 2336 2337 #undef __FUNCT__ 2338 #define __FUNCT__ "SNESGetFunction" 2339 /*@C 2340 SNESGetFunction - Returns the vector where the function is stored. 2341 2342 Not Collective, but Vec is parallel if SNES is parallel 2343 2344 Input Parameter: 2345 . snes - the SNES context 2346 2347 Output Parameter: 2348 + r - the function (or PETSC_NULL) 2349 . func - the function (or PETSC_NULL) 2350 - ctx - the function context (or PETSC_NULL) 2351 2352 Level: advanced 2353 2354 .keywords: SNES, nonlinear, get, function 2355 2356 .seealso: SNESSetFunction(), SNESGetSolution() 2357 @*/ 2358 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetFunction(SNES snes,Vec *r,PetscErrorCode (**func)(SNES,Vec,Vec,void*),void **ctx) 2359 { 2360 PetscFunctionBegin; 2361 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2362 if (r) *r = snes->vec_func; 2363 if (func) *func = snes->ops->computefunction; 2364 if (ctx) *ctx = snes->funP; 2365 PetscFunctionReturn(0); 2366 } 2367 2368 #undef __FUNCT__ 2369 #define __FUNCT__ "SNESSetOptionsPrefix" 2370 /*@C 2371 SNESSetOptionsPrefix - Sets the prefix used for searching for all 2372 SNES options in the database. 2373 2374 Collective on SNES 2375 2376 Input Parameter: 2377 + snes - the SNES context 2378 - prefix - the prefix to prepend to all option names 2379 2380 Notes: 2381 A hyphen (-) must NOT be given at the beginning of the prefix name. 2382 The first character of all runtime options is AUTOMATICALLY the hyphen. 2383 2384 Level: advanced 2385 2386 .keywords: SNES, set, options, prefix, database 2387 2388 .seealso: SNESSetFromOptions() 2389 @*/ 2390 PetscErrorCode PETSCSNES_DLLEXPORT SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2391 { 2392 PetscErrorCode ierr; 2393 2394 PetscFunctionBegin; 2395 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2396 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2397 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2398 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2399 PetscFunctionReturn(0); 2400 } 2401 2402 #undef __FUNCT__ 2403 #define __FUNCT__ "SNESAppendOptionsPrefix" 2404 /*@C 2405 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2406 SNES options in the database. 2407 2408 Collective on SNES 2409 2410 Input Parameters: 2411 + snes - the SNES context 2412 - prefix - the prefix to prepend to all option names 2413 2414 Notes: 2415 A hyphen (-) must NOT be given at the beginning of the prefix name. 2416 The first character of all runtime options is AUTOMATICALLY the hyphen. 2417 2418 Level: advanced 2419 2420 .keywords: SNES, append, options, prefix, database 2421 2422 .seealso: SNESGetOptionsPrefix() 2423 @*/ 2424 PetscErrorCode PETSCSNES_DLLEXPORT SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2425 { 2426 PetscErrorCode ierr; 2427 2428 PetscFunctionBegin; 2429 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2430 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2431 if (!snes->ksp) {ierr = SNESGetKSP(snes,&snes->ksp);CHKERRQ(ierr);} 2432 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2433 PetscFunctionReturn(0); 2434 } 2435 2436 #undef __FUNCT__ 2437 #define __FUNCT__ "SNESGetOptionsPrefix" 2438 /*@C 2439 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2440 SNES options in the database. 2441 2442 Not Collective 2443 2444 Input Parameter: 2445 . snes - the SNES context 2446 2447 Output Parameter: 2448 . prefix - pointer to the prefix string used 2449 2450 Notes: On the fortran side, the user should pass in a string 'prefix' of 2451 sufficient length to hold the prefix. 2452 2453 Level: advanced 2454 2455 .keywords: SNES, get, options, prefix, database 2456 2457 .seealso: SNESAppendOptionsPrefix() 2458 @*/ 2459 PetscErrorCode PETSCSNES_DLLEXPORT SNESGetOptionsPrefix(SNES snes,const char *prefix[]) 2460 { 2461 PetscErrorCode ierr; 2462 2463 PetscFunctionBegin; 2464 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2465 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2466 PetscFunctionReturn(0); 2467 } 2468 2469 2470 #undef __FUNCT__ 2471 #define __FUNCT__ "SNESRegister" 2472 /*@C 2473 SNESRegister - See SNESRegisterDynamic() 2474 2475 Level: advanced 2476 @*/ 2477 PetscErrorCode PETSCSNES_DLLEXPORT SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2478 { 2479 char fullname[PETSC_MAX_PATH_LEN]; 2480 PetscErrorCode ierr; 2481 2482 PetscFunctionBegin; 2483 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2484 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2485 PetscFunctionReturn(0); 2486 } 2487 2488 #undef __FUNCT__ 2489 #define __FUNCT__ "SNESTestLocalMin" 2490 PetscErrorCode PETSCSNES_DLLEXPORT SNESTestLocalMin(SNES snes) 2491 { 2492 PetscErrorCode ierr; 2493 PetscInt N,i,j; 2494 Vec u,uh,fh; 2495 PetscScalar value; 2496 PetscReal norm; 2497 2498 PetscFunctionBegin; 2499 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2500 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2501 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2502 2503 /* currently only works for sequential */ 2504 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2505 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2506 for (i=0; i<N; i++) { 2507 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2508 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2509 for (j=-10; j<11; j++) { 2510 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2511 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2512 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2513 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2514 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2515 value = -value; 2516 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2517 } 2518 } 2519 ierr = VecDestroy(uh);CHKERRQ(ierr); 2520 ierr = VecDestroy(fh);CHKERRQ(ierr); 2521 PetscFunctionReturn(0); 2522 } 2523 2524 #undef __FUNCT__ 2525 #define __FUNCT__ "SNESKSPSetUseEW" 2526 /*@ 2527 SNESKSPSetUseEW - Sets SNES use Eisenstat-Walker method for 2528 computing relative tolerance for linear solvers within an inexact 2529 Newton method. 2530 2531 Collective on SNES 2532 2533 Input Parameters: 2534 + snes - SNES context 2535 - flag - PETSC_TRUE or PETSC_FALSE 2536 2537 Notes: 2538 Currently, the default is to use a constant relative tolerance for 2539 the inner linear solvers. Alternatively, one can use the 2540 Eisenstat-Walker method, where the relative convergence tolerance 2541 is reset at each Newton iteration according progress of the nonlinear 2542 solver. 2543 2544 Level: advanced 2545 2546 Reference: 2547 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2548 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 2549 2550 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 2551 2552 .seealso: SNESKSPGetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 2553 @*/ 2554 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetUseEW(SNES snes,PetscTruth flag) 2555 { 2556 PetscFunctionBegin; 2557 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2558 snes->ksp_ewconv = flag; 2559 PetscFunctionReturn(0); 2560 } 2561 2562 #undef __FUNCT__ 2563 #define __FUNCT__ "SNESKSPGetUseEW" 2564 /*@ 2565 SNESKSPGetUseEW - Gets if SNES is using Eisenstat-Walker method 2566 for computing relative tolerance for linear solvers within an 2567 inexact Newton method. 2568 2569 Not Collective 2570 2571 Input Parameter: 2572 . snes - SNES context 2573 2574 Output Parameter: 2575 . flag - PETSC_TRUE or PETSC_FALSE 2576 2577 Notes: 2578 Currently, the default is to use a constant relative tolerance for 2579 the inner linear solvers. Alternatively, one can use the 2580 Eisenstat-Walker method, where the relative convergence tolerance 2581 is reset at each Newton iteration according progress of the nonlinear 2582 solver. 2583 2584 Level: advanced 2585 2586 Reference: 2587 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2588 inexact Newton method", SISC 17 (1), pp.16-32, 1996. 2589 2590 .keywords: SNES, KSP, Eisenstat, Walker, convergence, test, inexact, Newton 2591 2592 .seealso: SNESKSPSetUseEW(), SNESKSPGetParametersEW(), SNESKSPSetParametersEW() 2593 @*/ 2594 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetUseEW(SNES snes, PetscTruth *flag) 2595 { 2596 PetscFunctionBegin; 2597 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2598 PetscValidPointer(flag,2); 2599 *flag = snes->ksp_ewconv; 2600 PetscFunctionReturn(0); 2601 } 2602 2603 #undef __FUNCT__ 2604 #define __FUNCT__ "SNESKSPSetParametersEW" 2605 /*@ 2606 SNESKSPSetParametersEW - Sets parameters for Eisenstat-Walker 2607 convergence criteria for the linear solvers within an inexact 2608 Newton method. 2609 2610 Collective on SNES 2611 2612 Input Parameters: 2613 + snes - SNES context 2614 . version - version 1, 2 (default is 2) or 3 2615 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 2616 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 2617 . gamma - multiplicative factor for version 2 rtol computation 2618 (0 <= gamma2 <= 1) 2619 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 2620 . alpha2 - power for safeguard 2621 - threshold - threshold for imposing safeguard (0 < threshold < 1) 2622 2623 Note: 2624 Version 3 was contributed by Luis Chacon, June 2006. 2625 2626 Use PETSC_DEFAULT to retain the default for any of the parameters. 2627 2628 Level: advanced 2629 2630 Reference: 2631 S. C. Eisenstat and H. F. Walker, "Choosing the forcing terms in an 2632 inexact Newton method", Utah State University Math. Stat. Dept. Res. 2633 Report 6/94/75, June, 1994, to appear in SIAM J. Sci. Comput. 2634 2635 .keywords: SNES, KSP, Eisenstat, Walker, set, parameters 2636 2637 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPGetParametersEW() 2638 @*/ 2639 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPSetParametersEW(SNES snes,PetscInt version,PetscReal rtol_0,PetscReal rtol_max, 2640 PetscReal gamma,PetscReal alpha,PetscReal alpha2,PetscReal threshold) 2641 { 2642 SNESKSPEW *kctx; 2643 PetscFunctionBegin; 2644 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2645 kctx = (SNESKSPEW*)snes->kspconvctx; 2646 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 2647 2648 if (version != PETSC_DEFAULT) kctx->version = version; 2649 if (rtol_0 != PETSC_DEFAULT) kctx->rtol_0 = rtol_0; 2650 if (rtol_max != PETSC_DEFAULT) kctx->rtol_max = rtol_max; 2651 if (gamma != PETSC_DEFAULT) kctx->gamma = gamma; 2652 if (alpha != PETSC_DEFAULT) kctx->alpha = alpha; 2653 if (alpha2 != PETSC_DEFAULT) kctx->alpha2 = alpha2; 2654 if (threshold != PETSC_DEFAULT) kctx->threshold = threshold; 2655 2656 if (kctx->version < 1 || kctx->version > 3) { 2657 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 and 3 are supported: %D",kctx->version); 2658 } 2659 if (kctx->rtol_0 < 0.0 || kctx->rtol_0 >= 1.0) { 2660 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_0 < 1.0: %G",kctx->rtol_0); 2661 } 2662 if (kctx->rtol_max < 0.0 || kctx->rtol_max >= 1.0) { 2663 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= rtol_max (%G) < 1.0\n",kctx->rtol_max); 2664 } 2665 if (kctx->gamma < 0.0 || kctx->gamma > 1.0) { 2666 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 <= gamma (%G) <= 1.0\n",kctx->gamma); 2667 } 2668 if (kctx->alpha <= 1.0 || kctx->alpha > 2.0) { 2669 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"1.0 < alpha (%G) <= 2.0\n",kctx->alpha); 2670 } 2671 if (kctx->threshold <= 0.0 || kctx->threshold >= 1.0) { 2672 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"0.0 < threshold (%G) < 1.0\n",kctx->threshold); 2673 } 2674 PetscFunctionReturn(0); 2675 } 2676 2677 #undef __FUNCT__ 2678 #define __FUNCT__ "SNESKSPGetParametersEW" 2679 /*@ 2680 SNESKSPGetParametersEW - Gets parameters for Eisenstat-Walker 2681 convergence criteria for the linear solvers within an inexact 2682 Newton method. 2683 2684 Not Collective 2685 2686 Input Parameters: 2687 snes - SNES context 2688 2689 Output Parameters: 2690 + version - version 1, 2 (default is 2) or 3 2691 . rtol_0 - initial relative tolerance (0 <= rtol_0 < 1) 2692 . rtol_max - maximum relative tolerance (0 <= rtol_max < 1) 2693 . gamma - multiplicative factor for version 2 rtol computation 2694 (0 <= gamma2 <= 1) 2695 . alpha - power for version 2 rtol computation (1 < alpha <= 2) 2696 . alpha2 - power for safeguard 2697 - threshold - threshold for imposing safeguard (0 < threshold < 1) 2698 2699 Level: advanced 2700 2701 .keywords: SNES, KSP, Eisenstat, Walker, get, parameters 2702 2703 .seealso: SNESKSPSetUseEW(), SNESKSPGetUseEW(), SNESKSPSetParametersEW() 2704 @*/ 2705 PetscErrorCode PETSCSNES_DLLEXPORT SNESKSPGetParametersEW(SNES snes,PetscInt *version,PetscReal *rtol_0,PetscReal *rtol_max, 2706 PetscReal *gamma,PetscReal *alpha,PetscReal *alpha2,PetscReal *threshold) 2707 { 2708 SNESKSPEW *kctx; 2709 PetscFunctionBegin; 2710 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2711 kctx = (SNESKSPEW*)snes->kspconvctx; 2712 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context existing"); 2713 if(version) *version = kctx->version; 2714 if(rtol_0) *rtol_0 = kctx->rtol_0; 2715 if(rtol_max) *rtol_max = kctx->rtol_max; 2716 if(gamma) *gamma = kctx->gamma; 2717 if(alpha) *alpha = kctx->alpha; 2718 if(alpha2) *alpha2 = kctx->alpha2; 2719 if(threshold) *threshold = kctx->threshold; 2720 PetscFunctionReturn(0); 2721 } 2722 2723 #undef __FUNCT__ 2724 #define __FUNCT__ "SNESKSPEW_PreSolve" 2725 static PetscErrorCode SNESKSPEW_PreSolve(SNES snes, KSP ksp, Vec b, Vec x) 2726 { 2727 PetscErrorCode ierr; 2728 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 2729 PetscReal rtol=PETSC_DEFAULT,stol; 2730 2731 PetscFunctionBegin; 2732 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 2733 if (!snes->iter) { /* first time in, so use the original user rtol */ 2734 rtol = kctx->rtol_0; 2735 } else { 2736 if (kctx->version == 1) { 2737 rtol = (snes->norm - kctx->lresid_last)/kctx->norm_last; 2738 if (rtol < 0.0) rtol = -rtol; 2739 stol = pow(kctx->rtol_last,kctx->alpha2); 2740 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 2741 } else if (kctx->version == 2) { 2742 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 2743 stol = kctx->gamma * pow(kctx->rtol_last,kctx->alpha); 2744 if (stol > kctx->threshold) rtol = PetscMax(rtol,stol); 2745 } else if (kctx->version == 3) {/* contributed by Luis Chacon, June 2006. */ 2746 rtol = kctx->gamma * pow(snes->norm/kctx->norm_last,kctx->alpha); 2747 /* safeguard: avoid sharp decrease of rtol */ 2748 stol = kctx->gamma*pow(kctx->rtol_last,kctx->alpha); 2749 stol = PetscMax(rtol,stol); 2750 rtol = PetscMin(kctx->rtol_0,stol); 2751 /* safeguard: avoid oversolving */ 2752 stol = kctx->gamma*(snes->ttol)/snes->norm; 2753 stol = PetscMax(rtol,stol); 2754 rtol = PetscMin(kctx->rtol_0,stol); 2755 } else SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE,"Only versions 1, 2 or 3 are supported: %D",kctx->version); 2756 } 2757 /* safeguard: avoid rtol greater than one */ 2758 rtol = PetscMin(rtol,kctx->rtol_max); 2759 ierr = KSPSetTolerances(ksp,rtol,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);CHKERRQ(ierr); 2760 ierr = PetscInfo3(snes,"iter %D, Eisenstat-Walker (version %D) KSP rtol=%G\n",snes->iter,kctx->version,rtol);CHKERRQ(ierr); 2761 PetscFunctionReturn(0); 2762 } 2763 2764 #undef __FUNCT__ 2765 #define __FUNCT__ "SNESKSPEW_PostSolve" 2766 static PetscErrorCode SNESKSPEW_PostSolve(SNES snes, KSP ksp, Vec b, Vec x) 2767 { 2768 PetscErrorCode ierr; 2769 SNESKSPEW *kctx = (SNESKSPEW*)snes->kspconvctx; 2770 PCSide pcside; 2771 Vec lres; 2772 2773 PetscFunctionBegin; 2774 if (!kctx) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"No Eisenstat-Walker context exists"); 2775 ierr = KSPGetTolerances(ksp,&kctx->rtol_last,0,0,0);CHKERRQ(ierr); 2776 ierr = SNESGetFunctionNorm(snes,&kctx->norm_last);CHKERRQ(ierr); 2777 if (kctx->version == 1) { 2778 ierr = KSPGetPreconditionerSide(ksp,&pcside);CHKERRQ(ierr); 2779 if (pcside == PC_RIGHT) { /* XXX Should we also test KSP_UNPRECONDITIONED_NORM ? */ 2780 /* KSP residual is true linear residual */ 2781 ierr = KSPGetResidualNorm(ksp,&kctx->lresid_last);CHKERRQ(ierr); 2782 } else { 2783 /* KSP residual is preconditioned residual */ 2784 /* compute true linear residual norm */ 2785 ierr = VecDuplicate(b,&lres);CHKERRQ(ierr); 2786 ierr = MatMult(snes->jacobian,x,lres);CHKERRQ(ierr); 2787 ierr = VecAYPX(lres,-1.0,b);CHKERRQ(ierr); 2788 ierr = VecNorm(lres,NORM_2,&kctx->lresid_last);CHKERRQ(ierr); 2789 ierr = VecDestroy(lres);CHKERRQ(ierr); 2790 } 2791 } 2792 PetscFunctionReturn(0); 2793 } 2794 2795 #undef __FUNCT__ 2796 #define __FUNCT__ "SNES_KSPSolve" 2797 PetscErrorCode SNES_KSPSolve(SNES snes, KSP ksp, Vec b, Vec x) 2798 { 2799 PetscErrorCode ierr; 2800 2801 PetscFunctionBegin; 2802 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PreSolve(snes,ksp,b,x);CHKERRQ(ierr); } 2803 ierr = KSPSolve(ksp,b,x);CHKERRQ(ierr); 2804 if (snes->ksp_ewconv) { ierr = SNESKSPEW_PostSolve(snes,ksp,b,x);CHKERRQ(ierr); } 2805 PetscFunctionReturn(0); 2806 } 2807