1 2 #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 3 4 PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 5 PetscFList SNESList = PETSC_NULL; 6 7 /* Logging support */ 8 PetscCookie SNES_COOKIE = 0; 9 PetscEvent SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0; 10 11 #undef __FUNCT__ 12 #define __FUNCT__ "SNESView" 13 /*@C 14 SNESView - Prints the SNES data structure. 15 16 Collective on SNES 17 18 Input Parameters: 19 + SNES - the SNES context 20 - viewer - visualization context 21 22 Options Database Key: 23 . -snes_view - Calls SNESView() at end of SNESSolve() 24 25 Notes: 26 The available visualization contexts include 27 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 28 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 29 output where only the first processor opens 30 the file. All other processors send their 31 data to the first processor to print. 32 33 The user can open an alternative visualization context with 34 PetscViewerASCIIOpen() - output to a specified file. 35 36 Level: beginner 37 38 .keywords: SNES, view 39 40 .seealso: PetscViewerASCIIOpen() 41 @*/ 42 PetscErrorCode SNESView(SNES snes,PetscViewer viewer) 43 { 44 SNES_KSP_EW_ConvCtx *kctx; 45 PetscErrorCode ierr; 46 KSP ksp; 47 char *type; 48 PetscTruth iascii,isstring; 49 50 PetscFunctionBegin; 51 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 52 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm); 53 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 54 PetscCheckSameComm(snes,1,viewer,2); 55 56 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);CHKERRQ(ierr); 57 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 58 if (iascii) { 59 if (snes->prefix) { 60 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr); 61 } else { 62 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 63 } 64 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 65 if (type) { 66 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 67 } else { 68 ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 69 } 70 if (snes->view) { 71 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 72 ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 73 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 74 } 75 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%D, maximum function evaluations=%D\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 76 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n", 77 snes->rtol,snes->abstol,snes->xtol);CHKERRQ(ierr); 78 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%D\n",snes->linear_its);CHKERRQ(ierr); 79 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%D\n",snes->nfuncs);CHKERRQ(ierr); 80 if (snes->ksp_ewconv) { 81 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 82 if (kctx) { 83 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %D)\n",kctx->version);CHKERRQ(ierr); 84 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 85 ierr = PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 86 } 87 } 88 } else if (isstring) { 89 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 90 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 91 } 92 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 93 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 94 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 95 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 96 PetscFunctionReturn(0); 97 } 98 99 /* 100 We retain a list of functions that also take SNES command 101 line options. These are called at the end SNESSetFromOptions() 102 */ 103 #define MAXSETFROMOPTIONS 5 104 static PetscInt numberofsetfromoptions; 105 static PetscErrorCode (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 106 107 #undef __FUNCT__ 108 #define __FUNCT__ "SNESAddOptionsChecker" 109 /*@C 110 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 111 112 Not Collective 113 114 Input Parameter: 115 . snescheck - function that checks for options 116 117 Level: developer 118 119 .seealso: SNESSetFromOptions() 120 @*/ 121 PetscErrorCode SNESAddOptionsChecker(PetscErrorCode (*snescheck)(SNES)) 122 { 123 PetscFunctionBegin; 124 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 125 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %D allowed", MAXSETFROMOPTIONS); 126 } 127 othersetfromoptions[numberofsetfromoptions++] = snescheck; 128 PetscFunctionReturn(0); 129 } 130 131 #undef __FUNCT__ 132 #define __FUNCT__ "SNESSetFromOptions" 133 /*@ 134 SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 135 136 Collective on SNES 137 138 Input Parameter: 139 . snes - the SNES context 140 141 Options Database Keys: 142 + -snes_type <type> - ls, tr, umls, umtr, test 143 . -snes_stol - convergence tolerance in terms of the norm 144 of the change in the solution between steps 145 . -snes_atol <abstol> - absolute tolerance of residual norm 146 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 147 . -snes_max_it <max_it> - maximum number of iterations 148 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 149 . -snes_max_fail <max_fail> - maximum number of failures 150 . -snes_trtol <trtol> - trust region tolerance 151 . -snes_no_convergence_test - skip convergence test in nonlinear or minimization 152 solver; hence iterations will continue until max_it 153 or some other criterion is reached. Saves expense 154 of convergence test 155 . -snes_monitor - prints residual norm at each iteration 156 . -snes_vecmonitor - plots solution at each iteration 157 . -snes_vecmonitor_update - plots update to solution at each iteration 158 . -snes_xmonitor - plots residual norm at each iteration 159 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 160 . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 161 - -snes_print_converged_reason - print the reason for convergence/divergence after each solve 162 163 Options Database for Eisenstat-Walker method: 164 + -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence 165 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 166 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 167 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 168 . -snes_ksp_ew_gamma <gamma> - Sets gamma 169 . -snes_ksp_ew_alpha <alpha> - Sets alpha 170 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 171 - -snes_ksp_ew_threshold <threshold> - Sets threshold 172 173 Notes: 174 To see all options, run your program with the -help option or consult 175 the users manual. 176 177 Level: beginner 178 179 .keywords: SNES, nonlinear, set, options, database 180 181 .seealso: SNESSetOptionsPrefix() 182 @*/ 183 PetscErrorCode SNESSetFromOptions(SNES snes) 184 { 185 KSP ksp; 186 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 187 PetscTruth flg; 188 PetscErrorCode ierr; 189 PetscInt i; 190 const char *deft; 191 char type[256]; 192 193 PetscFunctionBegin; 194 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 195 196 ierr = PetscOptionsBegin(snes->comm,snes->prefix,"Nonlinear solver (SNES) options","SNES");CHKERRQ(ierr); 197 if (snes->type_name) { 198 deft = snes->type_name; 199 } else { 200 deft = SNESLS; 201 } 202 203 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 204 ierr = PetscOptionsList("-snes_type","Nonlinear solver method","SNESSetType",SNESList,deft,type,256,&flg);CHKERRQ(ierr); 205 if (flg) { 206 ierr = SNESSetType(snes,type);CHKERRQ(ierr); 207 } else if (!snes->type_name) { 208 ierr = SNESSetType(snes,deft);CHKERRQ(ierr); 209 } 210 ierr = PetscOptionsName("-snes_view","Print detailed information on solver used","SNESView",0);CHKERRQ(ierr); 211 212 ierr = PetscOptionsReal("-snes_stol","Stop if step length less then","SNESSetTolerances",snes->xtol,&snes->xtol,0);CHKERRQ(ierr); 213 ierr = PetscOptionsReal("-snes_atol","Stop if function norm less then","SNESSetTolerances",snes->abstol,&snes->abstol,0);CHKERRQ(ierr); 214 215 ierr = PetscOptionsReal("-snes_rtol","Stop if decrease in function norm less then","SNESSetTolerances",snes->rtol,&snes->rtol,0);CHKERRQ(ierr); 216 ierr = PetscOptionsInt("-snes_max_it","Maximum iterations","SNESSetTolerances",snes->max_its,&snes->max_its,PETSC_NULL);CHKERRQ(ierr); 217 ierr = PetscOptionsInt("-snes_max_funcs","Maximum function evaluations","SNESSetTolerances",snes->max_funcs,&snes->max_funcs,PETSC_NULL);CHKERRQ(ierr); 218 ierr = PetscOptionsInt("-snes_max_fail","Maximum failures","SNESSetTolerances",snes->maxFailures,&snes->maxFailures,PETSC_NULL);CHKERRQ(ierr); 219 ierr = PetscOptionsName("-snes_converged_reason","Print reason for converged or diverged","SNESSolve",&flg);CHKERRQ(ierr); 220 if (flg) { 221 snes->printreason = PETSC_TRUE; 222 } 223 224 ierr = PetscOptionsName("-snes_ksp_ew_conv","Use Eisentat-Walker linear system convergence test","SNES_KSP_SetParametersEW",&snes->ksp_ewconv);CHKERRQ(ierr); 225 226 ierr = PetscOptionsInt("-snes_ksp_ew_version","Version 1 or 2","SNES_KSP_SetParametersEW",kctx->version,&kctx->version,0);CHKERRQ(ierr); 227 ierr = PetscOptionsReal("-snes_ksp_ew_rtol0","0 <= rtol0 < 1","SNES_KSP_SetParametersEW",kctx->rtol_0,&kctx->rtol_0,0);CHKERRQ(ierr); 228 ierr = PetscOptionsReal("-snes_ksp_ew_rtolmax","0 <= rtolmax < 1","SNES_KSP_SetParametersEW",kctx->rtol_max,&kctx->rtol_max,0);CHKERRQ(ierr); 229 ierr = PetscOptionsReal("-snes_ksp_ew_gamma","0 <= gamma <= 1","SNES_KSP_SetParametersEW",kctx->gamma,&kctx->gamma,0);CHKERRQ(ierr); 230 ierr = PetscOptionsReal("-snes_ksp_ew_alpha","1 < alpha <= 2","SNES_KSP_SetParametersEW",kctx->alpha,&kctx->alpha,0);CHKERRQ(ierr); 231 ierr = PetscOptionsReal("-snes_ksp_ew_alpha2","alpha2","SNES_KSP_SetParametersEW",kctx->alpha2,&kctx->alpha2,0);CHKERRQ(ierr); 232 ierr = PetscOptionsReal("-snes_ksp_ew_threshold","0 < threshold < 1","SNES_KSP_SetParametersEW",kctx->threshold,&kctx->threshold,0);CHKERRQ(ierr); 233 234 ierr = PetscOptionsName("-snes_no_convergence_test","Don't test for convergence","None",&flg);CHKERRQ(ierr); 235 if (flg) {snes->converged = 0;} 236 ierr = PetscOptionsName("-snes_cancelmonitors","Remove all monitors","SNESClearMonitor",&flg);CHKERRQ(ierr); 237 if (flg) {ierr = SNESClearMonitor(snes);CHKERRQ(ierr);} 238 ierr = PetscOptionsName("-snes_monitor","Monitor norm of function","SNESDefaultMonitor",&flg);CHKERRQ(ierr); 239 if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultMonitor,0,0);CHKERRQ(ierr);} 240 ierr = PetscOptionsName("-snes_ratiomonitor","Monitor norm of function","SNESSetRatioMonitor",&flg);CHKERRQ(ierr); 241 if (flg) {ierr = SNESSetRatioMonitor(snes);CHKERRQ(ierr);} 242 ierr = PetscOptionsName("-snes_smonitor","Monitor norm of function (fewer digits)","SNESDefaultSMonitor",&flg);CHKERRQ(ierr); 243 if (flg) {ierr = SNESSetMonitor(snes,SNESDefaultSMonitor,0,0);CHKERRQ(ierr);} 244 ierr = PetscOptionsName("-snes_vecmonitor","Plot solution at each iteration","SNESVecViewMonitor",&flg);CHKERRQ(ierr); 245 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewMonitor,0,0);CHKERRQ(ierr);} 246 ierr = PetscOptionsName("-snes_vecmonitor_update","Plot correction at each iteration","SNESVecViewUpdateMonitor",&flg);CHKERRQ(ierr); 247 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewUpdateMonitor,0,0);CHKERRQ(ierr);} 248 ierr = PetscOptionsName("-snes_vecmonitor_residual","Plot residual at each iteration","SNESVecViewResidualMonitor",&flg);CHKERRQ(ierr); 249 if (flg) {ierr = SNESSetMonitor(snes,SNESVecViewResidualMonitor,0,0);CHKERRQ(ierr);} 250 ierr = PetscOptionsName("-snes_xmonitor","Plot function norm at each iteration","SNESLGMonitor",&flg);CHKERRQ(ierr); 251 if (flg) {ierr = SNESSetMonitor(snes,SNESLGMonitor,PETSC_NULL,PETSC_NULL);CHKERRQ(ierr);} 252 253 ierr = PetscOptionsName("-snes_fd","Use finite differences (slow) to compute Jacobian","SNESDefaultComputeJacobian",&flg);CHKERRQ(ierr); 254 if (flg) { 255 ierr = SNESSetJacobian(snes,snes->jacobian,snes->jacobian_pre,SNESDefaultComputeJacobian,snes->funP);CHKERRQ(ierr); 256 PetscLogInfo(snes,"SNESSetFromOptions: Setting default finite difference Jacobian matrix\n"); 257 } 258 259 for(i = 0; i < numberofsetfromoptions; i++) { 260 ierr = (*othersetfromoptions[i])(snes);CHKERRQ(ierr); 261 } 262 263 if (snes->setfromoptions) { 264 ierr = (*snes->setfromoptions)(snes);CHKERRQ(ierr); 265 } 266 267 ierr = PetscOptionsEnd();CHKERRQ(ierr); 268 269 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 270 ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 271 272 PetscFunctionReturn(0); 273 } 274 275 276 #undef __FUNCT__ 277 #define __FUNCT__ "SNESSetApplicationContext" 278 /*@ 279 SNESSetApplicationContext - Sets the optional user-defined context for 280 the nonlinear solvers. 281 282 Collective on SNES 283 284 Input Parameters: 285 + snes - the SNES context 286 - usrP - optional user context 287 288 Level: intermediate 289 290 .keywords: SNES, nonlinear, set, application, context 291 292 .seealso: SNESGetApplicationContext() 293 @*/ 294 PetscErrorCode SNESSetApplicationContext(SNES snes,void *usrP) 295 { 296 PetscFunctionBegin; 297 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 298 snes->user = usrP; 299 PetscFunctionReturn(0); 300 } 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "SNESGetApplicationContext" 304 /*@C 305 SNESGetApplicationContext - Gets the user-defined context for the 306 nonlinear solvers. 307 308 Not Collective 309 310 Input Parameter: 311 . snes - SNES context 312 313 Output Parameter: 314 . usrP - user context 315 316 Level: intermediate 317 318 .keywords: SNES, nonlinear, get, application, context 319 320 .seealso: SNESSetApplicationContext() 321 @*/ 322 PetscErrorCode SNESGetApplicationContext(SNES snes,void **usrP) 323 { 324 PetscFunctionBegin; 325 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 326 *usrP = snes->user; 327 PetscFunctionReturn(0); 328 } 329 330 #undef __FUNCT__ 331 #define __FUNCT__ "SNESGetIterationNumber" 332 /*@ 333 SNESGetIterationNumber - Gets the number of nonlinear iterations completed 334 at this time. 335 336 Not Collective 337 338 Input Parameter: 339 . snes - SNES context 340 341 Output Parameter: 342 . iter - iteration number 343 344 Notes: 345 For example, during the computation of iteration 2 this would return 1. 346 347 This is useful for using lagged Jacobians (where one does not recompute the 348 Jacobian at each SNES iteration). For example, the code 349 .vb 350 ierr = SNESGetIterationNumber(snes,&it); 351 if (!(it % 2)) { 352 [compute Jacobian here] 353 } 354 .ve 355 can be used in your ComputeJacobian() function to cause the Jacobian to be 356 recomputed every second SNES iteration. 357 358 Level: intermediate 359 360 .keywords: SNES, nonlinear, get, iteration, number 361 @*/ 362 PetscErrorCode SNESGetIterationNumber(SNES snes,PetscInt* iter) 363 { 364 PetscFunctionBegin; 365 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 366 PetscValidIntPointer(iter,2); 367 *iter = snes->iter; 368 PetscFunctionReturn(0); 369 } 370 371 #undef __FUNCT__ 372 #define __FUNCT__ "SNESGetFunctionNorm" 373 /*@ 374 SNESGetFunctionNorm - Gets the norm of the current function that was set 375 with SNESSSetFunction(). 376 377 Collective on SNES 378 379 Input Parameter: 380 . snes - SNES context 381 382 Output Parameter: 383 . fnorm - 2-norm of function 384 385 Level: intermediate 386 387 .keywords: SNES, nonlinear, get, function, norm 388 389 .seealso: SNESGetFunction() 390 @*/ 391 PetscErrorCode SNESGetFunctionNorm(SNES snes,PetscScalar *fnorm) 392 { 393 PetscFunctionBegin; 394 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 395 PetscValidScalarPointer(fnorm,2); 396 *fnorm = snes->norm; 397 PetscFunctionReturn(0); 398 } 399 400 #undef __FUNCT__ 401 #define __FUNCT__ "SNESGetNumberUnsuccessfulSteps" 402 /*@ 403 SNESGetNumberUnsuccessfulSteps - Gets the number of unsuccessful steps 404 attempted by the nonlinear solver. 405 406 Not Collective 407 408 Input Parameter: 409 . snes - SNES context 410 411 Output Parameter: 412 . nfails - number of unsuccessful steps attempted 413 414 Notes: 415 This counter is reset to zero for each successive call to SNESSolve(). 416 417 Level: intermediate 418 419 .keywords: SNES, nonlinear, get, number, unsuccessful, steps 420 @*/ 421 PetscErrorCode SNESGetNumberUnsuccessfulSteps(SNES snes,PetscInt* nfails) 422 { 423 PetscFunctionBegin; 424 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 425 PetscValidIntPointer(nfails,2); 426 *nfails = snes->numFailures; 427 PetscFunctionReturn(0); 428 } 429 430 #undef __FUNCT__ 431 #define __FUNCT__ "SNESSetMaximumUnsuccessfulSteps" 432 /*@ 433 SNESSetMaximumUnsuccessfulSteps - Sets the maximum number of unsuccessful steps 434 attempted by the nonlinear solver before it gives up. 435 436 Not Collective 437 438 Input Parameters: 439 + snes - SNES context 440 - maxFails - maximum of unsuccessful steps 441 442 Level: intermediate 443 444 .keywords: SNES, nonlinear, set, maximum, unsuccessful, steps 445 @*/ 446 PetscErrorCode SNESSetMaximumUnsuccessfulSteps(SNES snes, PetscInt maxFails) 447 { 448 PetscFunctionBegin; 449 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 450 snes->maxFailures = maxFails; 451 PetscFunctionReturn(0); 452 } 453 454 #undef __FUNCT__ 455 #define __FUNCT__ "SNESGetMaximumUnsuccessfulSteps" 456 /*@ 457 SNESGetMaximumUnsuccessfulSteps - Gets the maximum number of unsuccessful steps 458 attempted by the nonlinear solver before it gives up. 459 460 Not Collective 461 462 Input Parameter: 463 . snes - SNES context 464 465 Output Parameter: 466 . maxFails - maximum of unsuccessful steps 467 468 Level: intermediate 469 470 .keywords: SNES, nonlinear, get, maximum, unsuccessful, steps 471 @*/ 472 PetscErrorCode SNESGetMaximumUnsuccessfulSteps(SNES snes, PetscInt *maxFails) 473 { 474 PetscFunctionBegin; 475 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 476 PetscValidIntPointer(maxFails,2); 477 *maxFails = snes->maxFailures; 478 PetscFunctionReturn(0); 479 } 480 481 #undef __FUNCT__ 482 #define __FUNCT__ "SNESGetNumberLinearIterations" 483 /*@ 484 SNESGetNumberLinearIterations - Gets the total number of linear iterations 485 used by the nonlinear solver. 486 487 Not Collective 488 489 Input Parameter: 490 . snes - SNES context 491 492 Output Parameter: 493 . lits - number of linear iterations 494 495 Notes: 496 This counter is reset to zero for each successive call to SNESSolve(). 497 498 Level: intermediate 499 500 .keywords: SNES, nonlinear, get, number, linear, iterations 501 @*/ 502 PetscErrorCode SNESGetNumberLinearIterations(SNES snes,PetscInt* lits) 503 { 504 PetscFunctionBegin; 505 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 506 PetscValidIntPointer(lits,2); 507 *lits = snes->linear_its; 508 PetscFunctionReturn(0); 509 } 510 511 #undef __FUNCT__ 512 #define __FUNCT__ "SNESGetKSP" 513 /*@C 514 SNESGetKSP - Returns the KSP context for a SNES solver. 515 516 Not Collective, but if SNES object is parallel, then KSP object is parallel 517 518 Input Parameter: 519 . snes - the SNES context 520 521 Output Parameter: 522 . ksp - the KSP context 523 524 Notes: 525 The user can then directly manipulate the KSP context to set various 526 options, etc. Likewise, the user can then extract and manipulate the 527 KSP and PC contexts as well. 528 529 Level: beginner 530 531 .keywords: SNES, nonlinear, get, KSP, context 532 533 .seealso: KSPGetPC() 534 @*/ 535 PetscErrorCode SNESGetKSP(SNES snes,KSP *ksp) 536 { 537 PetscFunctionBegin; 538 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 539 PetscValidPointer(ksp,2); 540 *ksp = snes->ksp; 541 PetscFunctionReturn(0); 542 } 543 544 #undef __FUNCT__ 545 #define __FUNCT__ "SNESPublish_Petsc" 546 static PetscErrorCode SNESPublish_Petsc(PetscObject obj) 547 { 548 PetscFunctionBegin; 549 PetscFunctionReturn(0); 550 } 551 552 /* -----------------------------------------------------------*/ 553 #undef __FUNCT__ 554 #define __FUNCT__ "SNESCreate" 555 /*@C 556 SNESCreate - Creates a nonlinear solver context. 557 558 Collective on MPI_Comm 559 560 Input Parameters: 561 + comm - MPI communicator 562 563 Output Parameter: 564 . outsnes - the new SNES context 565 566 Options Database Keys: 567 + -snes_mf - Activates default matrix-free Jacobian-vector products, 568 and no preconditioning matrix 569 . -snes_mf_operator - Activates default matrix-free Jacobian-vector 570 products, and a user-provided preconditioning matrix 571 as set by SNESSetJacobian() 572 - -snes_fd - Uses (slow!) finite differences to compute Jacobian 573 574 Level: beginner 575 576 .keywords: SNES, nonlinear, create, context 577 578 .seealso: SNESSolve(), SNESDestroy(), SNES 579 @*/ 580 PetscErrorCode SNESCreate(MPI_Comm comm,SNES *outsnes) 581 { 582 PetscErrorCode ierr; 583 SNES snes; 584 SNES_KSP_EW_ConvCtx *kctx; 585 586 PetscFunctionBegin; 587 PetscValidPointer(outsnes,1); 588 *outsnes = PETSC_NULL; 589 #ifndef PETSC_USE_DYNAMIC_LIBRARIES 590 ierr = SNESInitializePackage(PETSC_NULL);CHKERRQ(ierr); 591 #endif 592 593 PetscHeaderCreate(snes,_p_SNES,PetscInt,SNES_COOKIE,0,"SNES",comm,SNESDestroy,SNESView); 594 PetscLogObjectCreate(snes); 595 snes->bops->publish = SNESPublish_Petsc; 596 snes->max_its = 50; 597 snes->max_funcs = 10000; 598 snes->norm = 0.0; 599 snes->rtol = 1.e-8; 600 snes->ttol = 0.0; 601 snes->abstol = 1.e-50; 602 snes->xtol = 1.e-8; 603 snes->deltatol = 1.e-12; 604 snes->nfuncs = 0; 605 snes->numFailures = 0; 606 snes->maxFailures = 1; 607 snes->linear_its = 0; 608 snes->numbermonitors = 0; 609 snes->data = 0; 610 snes->view = 0; 611 snes->setupcalled = 0; 612 snes->ksp_ewconv = PETSC_FALSE; 613 snes->vwork = 0; 614 snes->nwork = 0; 615 snes->conv_hist_len = 0; 616 snes->conv_hist_max = 0; 617 snes->conv_hist = PETSC_NULL; 618 snes->conv_hist_its = PETSC_NULL; 619 snes->conv_hist_reset = PETSC_TRUE; 620 snes->reason = SNES_CONVERGED_ITERATING; 621 622 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 623 ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 624 PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 625 snes->kspconvctx = (void*)kctx; 626 kctx->version = 2; 627 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 628 this was too large for some test cases */ 629 kctx->rtol_last = 0; 630 kctx->rtol_max = .9; 631 kctx->gamma = 1.0; 632 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 633 kctx->alpha = kctx->alpha2; 634 kctx->threshold = .1; 635 kctx->lresid_last = 0; 636 kctx->norm_last = 0; 637 638 ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr); 639 PetscLogObjectParent(snes,snes->ksp) 640 641 *outsnes = snes; 642 ierr = PetscPublishAll(snes);CHKERRQ(ierr); 643 PetscFunctionReturn(0); 644 } 645 646 #undef __FUNCT__ 647 #define __FUNCT__ "SNESSetFunction" 648 /*@C 649 SNESSetFunction - Sets the function evaluation routine and function 650 vector for use by the SNES routines in solving systems of nonlinear 651 equations. 652 653 Collective on SNES 654 655 Input Parameters: 656 + snes - the SNES context 657 . func - function evaluation routine 658 . r - vector to store function value 659 - ctx - [optional] user-defined context for private data for the 660 function evaluation routine (may be PETSC_NULL) 661 662 Calling sequence of func: 663 $ func (SNES snes,Vec x,Vec f,void *ctx); 664 665 . f - function vector 666 - ctx - optional user-defined function context 667 668 Notes: 669 The Newton-like methods typically solve linear systems of the form 670 $ f'(x) x = -f(x), 671 where f'(x) denotes the Jacobian matrix and f(x) is the function. 672 673 Level: beginner 674 675 .keywords: SNES, nonlinear, set, function 676 677 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 678 @*/ 679 PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 680 { 681 PetscFunctionBegin; 682 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 683 PetscValidHeaderSpecific(r,VEC_COOKIE,2); 684 PetscCheckSameComm(snes,1,r,2); 685 686 snes->computefunction = func; 687 snes->vec_func = snes->vec_func_always = r; 688 snes->funP = ctx; 689 PetscFunctionReturn(0); 690 } 691 692 /* --------------------------------------------------------------- */ 693 #undef __FUNCT__ 694 #define __FUNCT__ "SNESSetRhs" 695 /*@C 696 SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set 697 it assumes a zero right hand side. 698 699 Collective on SNES 700 701 Input Parameters: 702 + snes - the SNES context 703 - rhs - the right hand side vector or PETSC_NULL for a zero right hand side 704 705 Level: intermediate 706 707 .keywords: SNES, nonlinear, set, function, right hand side 708 709 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 710 @*/ 711 PetscErrorCode SNESSetRhs(SNES snes,Vec rhs) 712 { 713 PetscErrorCode ierr; 714 715 PetscFunctionBegin; 716 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 717 if (rhs) { 718 PetscValidHeaderSpecific(rhs,VEC_COOKIE,2); 719 PetscCheckSameComm(snes,1,rhs,2); 720 ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr); 721 } 722 if (snes->afine) { 723 ierr = VecDestroy(snes->afine);CHKERRQ(ierr); 724 } 725 snes->afine = rhs; 726 PetscFunctionReturn(0); 727 } 728 729 #undef __FUNCT__ 730 #define __FUNCT__ "SNESComputeFunction" 731 /*@ 732 SNESComputeFunction - Calls the function that has been set with 733 SNESSetFunction(). 734 735 Collective on SNES 736 737 Input Parameters: 738 + snes - the SNES context 739 - x - input vector 740 741 Output Parameter: 742 . y - function vector, as set by SNESSetFunction() 743 744 Notes: 745 SNESComputeFunction() is typically used within nonlinear solvers 746 implementations, so most users would not generally call this routine 747 themselves. 748 749 Level: developer 750 751 .keywords: SNES, nonlinear, compute, function 752 753 .seealso: SNESSetFunction(), SNESGetFunction() 754 @*/ 755 PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y) 756 { 757 PetscErrorCode ierr; 758 759 PetscFunctionBegin; 760 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 761 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 762 PetscValidHeaderSpecific(y,VEC_COOKIE,3); 763 PetscCheckSameComm(snes,1,x,2); 764 PetscCheckSameComm(snes,1,y,3); 765 766 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 767 PetscStackPush("SNES user function"); 768 ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 769 PetscStackPop; 770 if (snes->afine) { 771 PetscScalar mone = -1.0; 772 ierr = VecAXPY(&mone,snes->afine,y);CHKERRQ(ierr); 773 } 774 snes->nfuncs++; 775 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 776 PetscFunctionReturn(0); 777 } 778 779 #undef __FUNCT__ 780 #define __FUNCT__ "SNESComputeJacobian" 781 /*@ 782 SNESComputeJacobian - Computes the Jacobian matrix that has been 783 set with SNESSetJacobian(). 784 785 Collective on SNES and Mat 786 787 Input Parameters: 788 + snes - the SNES context 789 - x - input vector 790 791 Output Parameters: 792 + A - Jacobian matrix 793 . B - optional preconditioning matrix 794 - flag - flag indicating matrix structure 795 796 Notes: 797 Most users should not need to explicitly call this routine, as it 798 is used internally within the nonlinear solvers. 799 800 See KSPSetOperators() for important information about setting the 801 flag parameter. 802 803 Level: developer 804 805 .keywords: SNES, compute, Jacobian, matrix 806 807 .seealso: SNESSetJacobian(), KSPSetOperators() 808 @*/ 809 PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 810 { 811 PetscErrorCode ierr; 812 813 PetscFunctionBegin; 814 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 815 PetscValidHeaderSpecific(X,VEC_COOKIE,2); 816 PetscValidPointer(flg,5); 817 PetscCheckSameComm(snes,1,X,2); 818 if (!snes->computejacobian) PetscFunctionReturn(0); 819 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 820 *flg = DIFFERENT_NONZERO_PATTERN; 821 PetscStackPush("SNES user Jacobian function"); 822 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 823 PetscStackPop; 824 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 825 /* make sure user returned a correct Jacobian and preconditioner */ 826 PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 827 PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 828 PetscFunctionReturn(0); 829 } 830 831 #undef __FUNCT__ 832 #define __FUNCT__ "SNESSetJacobian" 833 /*@C 834 SNESSetJacobian - Sets the function to compute Jacobian as well as the 835 location to store the matrix. 836 837 Collective on SNES and Mat 838 839 Input Parameters: 840 + snes - the SNES context 841 . A - Jacobian matrix 842 . B - preconditioner matrix (usually same as the Jacobian) 843 . func - Jacobian evaluation routine 844 - ctx - [optional] user-defined context for private data for the 845 Jacobian evaluation routine (may be PETSC_NULL) 846 847 Calling sequence of func: 848 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 849 850 + x - input vector 851 . A - Jacobian matrix 852 . B - preconditioner matrix, usually the same as A 853 . flag - flag indicating information about the preconditioner matrix 854 structure (same as flag in KSPSetOperators()) 855 - ctx - [optional] user-defined Jacobian context 856 857 Notes: 858 See KSPSetOperators() for important information about setting the flag 859 output parameter in the routine func(). Be sure to read this information! 860 861 The routine func() takes Mat * as the matrix arguments rather than Mat. 862 This allows the Jacobian evaluation routine to replace A and/or B with a 863 completely new new matrix structure (not just different matrix elements) 864 when appropriate, for instance, if the nonzero structure is changing 865 throughout the global iterations. 866 867 Level: beginner 868 869 .keywords: SNES, nonlinear, set, Jacobian, matrix 870 871 .seealso: KSPSetOperators(), SNESSetFunction(), , MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor() 872 @*/ 873 PetscErrorCode SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 874 { 875 PetscErrorCode ierr; 876 877 PetscFunctionBegin; 878 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 879 if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 880 if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 881 if (A) PetscCheckSameComm(snes,1,A,2); 882 if (B) PetscCheckSameComm(snes,1,B,2); 883 if (func) snes->computejacobian = func; 884 if (ctx) snes->jacP = ctx; 885 if (A) { 886 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 887 snes->jacobian = A; 888 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 889 } 890 if (B) { 891 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 892 snes->jacobian_pre = B; 893 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 894 } 895 ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 896 PetscFunctionReturn(0); 897 } 898 899 #undef __FUNCT__ 900 #define __FUNCT__ "SNESGetJacobian" 901 /*@C 902 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 903 provided context for evaluating the Jacobian. 904 905 Not Collective, but Mat object will be parallel if SNES object is 906 907 Input Parameter: 908 . snes - the nonlinear solver context 909 910 Output Parameters: 911 + A - location to stash Jacobian matrix (or PETSC_NULL) 912 . B - location to stash preconditioner matrix (or PETSC_NULL) 913 . ctx - location to stash Jacobian ctx (or PETSC_NULL) 914 - func - location to put Jacobian function (or PETSC_NULL) 915 916 Level: advanced 917 918 .seealso: SNESSetJacobian(), SNESComputeJacobian() 919 @*/ 920 PetscErrorCode SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 921 { 922 PetscFunctionBegin; 923 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 924 if (A) *A = snes->jacobian; 925 if (B) *B = snes->jacobian_pre; 926 if (ctx) *ctx = snes->jacP; 927 if (func) *func = snes->computejacobian; 928 PetscFunctionReturn(0); 929 } 930 931 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 932 EXTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 933 934 #undef __FUNCT__ 935 #define __FUNCT__ "SNESSetUp" 936 /*@ 937 SNESSetUp - Sets up the internal data structures for the later use 938 of a nonlinear solver. 939 940 Collective on SNES 941 942 Input Parameters: 943 + snes - the SNES context 944 - x - the solution vector 945 946 Notes: 947 For basic use of the SNES solvers the user need not explicitly call 948 SNESSetUp(), since these actions will automatically occur during 949 the call to SNESSolve(). However, if one wishes to control this 950 phase separately, SNESSetUp() should be called after SNESCreate() 951 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 952 953 Level: advanced 954 955 .keywords: SNES, nonlinear, setup 956 957 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 958 @*/ 959 PetscErrorCode SNESSetUp(SNES snes,Vec x) 960 { 961 PetscErrorCode ierr; 962 PetscTruth flg, iseqtr; 963 964 PetscFunctionBegin; 965 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 966 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 967 PetscCheckSameComm(snes,1,x,2); 968 snes->vec_sol = snes->vec_sol_always = x; 969 970 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 971 /* 972 This version replaces the user provided Jacobian matrix with a 973 matrix-free version but still employs the user-provided preconditioner matrix 974 */ 975 if (flg) { 976 Mat J; 977 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 978 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 979 PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n"); 980 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 981 ierr = MatDestroy(J);CHKERRQ(ierr); 982 } 983 984 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) 985 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 986 if (flg) { 987 Mat J; 988 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 989 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 990 ierr = MatDestroy(J);CHKERRQ(ierr); 991 } 992 #endif 993 994 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 995 /* 996 This version replaces both the user-provided Jacobian and the user- 997 provided preconditioner matrix with the default matrix free version. 998 */ 999 if (flg) { 1000 Mat J; 1001 KSP ksp; 1002 PC pc; 1003 1004 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1005 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1006 PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"); 1007 ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 1008 ierr = MatDestroy(J);CHKERRQ(ierr); 1009 1010 /* force no preconditioner */ 1011 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1012 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1013 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1014 if (!flg) { 1015 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1016 } 1017 } 1018 1019 if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1020 if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1021 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1022 if (snes->vec_func == snes->vec_sol) { 1023 SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1024 } 1025 1026 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1027 ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 1028 if (snes->ksp_ewconv && !iseqtr) { 1029 KSP ksp; 1030 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1031 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1032 } 1033 1034 if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 1035 snes->setupcalled = 1; 1036 PetscFunctionReturn(0); 1037 } 1038 1039 #undef __FUNCT__ 1040 #define __FUNCT__ "SNESDestroy" 1041 /*@C 1042 SNESDestroy - Destroys the nonlinear solver context that was created 1043 with SNESCreate(). 1044 1045 Collective on SNES 1046 1047 Input Parameter: 1048 . snes - the SNES context 1049 1050 Level: beginner 1051 1052 .keywords: SNES, nonlinear, destroy 1053 1054 .seealso: SNESCreate(), SNESSolve() 1055 @*/ 1056 PetscErrorCode SNESDestroy(SNES snes) 1057 { 1058 PetscInt i; 1059 PetscErrorCode ierr; 1060 1061 PetscFunctionBegin; 1062 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1063 if (--snes->refct > 0) PetscFunctionReturn(0); 1064 1065 /* if memory was published with AMS then destroy it */ 1066 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1067 1068 if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1069 if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 1070 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1071 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1072 if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 1073 ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1074 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1075 for (i=0; i<snes->numbermonitors; i++) { 1076 if (snes->monitordestroy[i]) { 1077 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1078 } 1079 } 1080 PetscLogObjectDestroy((PetscObject)snes); 1081 PetscHeaderDestroy((PetscObject)snes); 1082 PetscFunctionReturn(0); 1083 } 1084 1085 /* ----------- Routines to set solver parameters ---------- */ 1086 1087 #undef __FUNCT__ 1088 #define __FUNCT__ "SNESSetTolerances" 1089 /*@ 1090 SNESSetTolerances - Sets various parameters used in convergence tests. 1091 1092 Collective on SNES 1093 1094 Input Parameters: 1095 + snes - the SNES context 1096 . abstol - absolute convergence tolerance 1097 . rtol - relative convergence tolerance 1098 . stol - convergence tolerance in terms of the norm 1099 of the change in the solution between steps 1100 . maxit - maximum number of iterations 1101 - maxf - maximum number of function evaluations 1102 1103 Options Database Keys: 1104 + -snes_atol <abstol> - Sets abstol 1105 . -snes_rtol <rtol> - Sets rtol 1106 . -snes_stol <stol> - Sets stol 1107 . -snes_max_it <maxit> - Sets maxit 1108 - -snes_max_funcs <maxf> - Sets maxf 1109 1110 Notes: 1111 The default maximum number of iterations is 50. 1112 The default maximum number of function evaluations is 1000. 1113 1114 Level: intermediate 1115 1116 .keywords: SNES, nonlinear, set, convergence, tolerances 1117 1118 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 1119 @*/ 1120 PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1121 { 1122 PetscFunctionBegin; 1123 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1124 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1125 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1126 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1127 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1128 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1129 PetscFunctionReturn(0); 1130 } 1131 1132 #undef __FUNCT__ 1133 #define __FUNCT__ "SNESGetTolerances" 1134 /*@ 1135 SNESGetTolerances - Gets various parameters used in convergence tests. 1136 1137 Not Collective 1138 1139 Input Parameters: 1140 + snes - the SNES context 1141 . abstol - absolute convergence tolerance 1142 . rtol - relative convergence tolerance 1143 . stol - convergence tolerance in terms of the norm 1144 of the change in the solution between steps 1145 . maxit - maximum number of iterations 1146 - maxf - maximum number of function evaluations 1147 1148 Notes: 1149 The user can specify PETSC_NULL for any parameter that is not needed. 1150 1151 Level: intermediate 1152 1153 .keywords: SNES, nonlinear, get, convergence, tolerances 1154 1155 .seealso: SNESSetTolerances() 1156 @*/ 1157 PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1158 { 1159 PetscFunctionBegin; 1160 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1161 if (abstol) *abstol = snes->abstol; 1162 if (rtol) *rtol = snes->rtol; 1163 if (stol) *stol = snes->xtol; 1164 if (maxit) *maxit = snes->max_its; 1165 if (maxf) *maxf = snes->max_funcs; 1166 PetscFunctionReturn(0); 1167 } 1168 1169 #undef __FUNCT__ 1170 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1171 /*@ 1172 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1173 1174 Collective on SNES 1175 1176 Input Parameters: 1177 + snes - the SNES context 1178 - tol - tolerance 1179 1180 Options Database Key: 1181 . -snes_trtol <tol> - Sets tol 1182 1183 Level: intermediate 1184 1185 .keywords: SNES, nonlinear, set, trust region, tolerance 1186 1187 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 1188 @*/ 1189 PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1190 { 1191 PetscFunctionBegin; 1192 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1193 snes->deltatol = tol; 1194 PetscFunctionReturn(0); 1195 } 1196 1197 /* 1198 Duplicate the lg monitors for SNES from KSP; for some reason with 1199 dynamic libraries things don't work under Sun4 if we just use 1200 macros instead of functions 1201 */ 1202 #undef __FUNCT__ 1203 #define __FUNCT__ "SNESLGMonitor" 1204 PetscErrorCode SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1205 { 1206 PetscErrorCode ierr; 1207 1208 PetscFunctionBegin; 1209 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1210 ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1211 PetscFunctionReturn(0); 1212 } 1213 1214 #undef __FUNCT__ 1215 #define __FUNCT__ "SNESLGMonitorCreate" 1216 PetscErrorCode SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1217 { 1218 PetscErrorCode ierr; 1219 1220 PetscFunctionBegin; 1221 ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1222 PetscFunctionReturn(0); 1223 } 1224 1225 #undef __FUNCT__ 1226 #define __FUNCT__ "SNESLGMonitorDestroy" 1227 PetscErrorCode SNESLGMonitorDestroy(PetscDrawLG draw) 1228 { 1229 PetscErrorCode ierr; 1230 1231 PetscFunctionBegin; 1232 ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1233 PetscFunctionReturn(0); 1234 } 1235 1236 /* ------------ Routines to set performance monitoring options ----------- */ 1237 1238 #undef __FUNCT__ 1239 #define __FUNCT__ "SNESSetMonitor" 1240 /*@C 1241 SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 1242 iteration of the nonlinear solver to display the iteration's 1243 progress. 1244 1245 Collective on SNES 1246 1247 Input Parameters: 1248 + snes - the SNES context 1249 . func - monitoring routine 1250 . mctx - [optional] user-defined context for private data for the 1251 monitor routine (use PETSC_NULL if no context is desitre) 1252 - monitordestroy - [optional] routine that frees monitor context 1253 (may be PETSC_NULL) 1254 1255 Calling sequence of func: 1256 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1257 1258 + snes - the SNES context 1259 . its - iteration number 1260 . norm - 2-norm function value (may be estimated) 1261 - mctx - [optional] monitoring context 1262 1263 Options Database Keys: 1264 + -snes_monitor - sets SNESDefaultMonitor() 1265 . -snes_xmonitor - sets line graph monitor, 1266 uses SNESLGMonitorCreate() 1267 _ -snes_cancelmonitors - cancels all monitors that have 1268 been hardwired into a code by 1269 calls to SNESSetMonitor(), but 1270 does not cancel those set via 1271 the options database. 1272 1273 Notes: 1274 Several different monitoring routines may be set by calling 1275 SNESSetMonitor() multiple times; all will be called in the 1276 order in which they were set. 1277 1278 Level: intermediate 1279 1280 .keywords: SNES, nonlinear, set, monitor 1281 1282 .seealso: SNESDefaultMonitor(), SNESClearMonitor() 1283 @*/ 1284 PetscErrorCode SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 1285 { 1286 PetscFunctionBegin; 1287 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1288 if (snes->numbermonitors >= MAXSNESMONITORS) { 1289 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1290 } 1291 snes->monitor[snes->numbermonitors] = func; 1292 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1293 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1294 PetscFunctionReturn(0); 1295 } 1296 1297 #undef __FUNCT__ 1298 #define __FUNCT__ "SNESClearMonitor" 1299 /*@C 1300 SNESClearMonitor - Clears all the monitor functions for a SNES object. 1301 1302 Collective on SNES 1303 1304 Input Parameters: 1305 . snes - the SNES context 1306 1307 Options Database Key: 1308 . -snes_cancelmonitors - cancels all monitors that have been hardwired 1309 into a code by calls to SNESSetMonitor(), but does not cancel those 1310 set via the options database 1311 1312 Notes: 1313 There is no way to clear one specific monitor from a SNES object. 1314 1315 Level: intermediate 1316 1317 .keywords: SNES, nonlinear, set, monitor 1318 1319 .seealso: SNESDefaultMonitor(), SNESSetMonitor() 1320 @*/ 1321 PetscErrorCode SNESClearMonitor(SNES snes) 1322 { 1323 PetscFunctionBegin; 1324 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1325 snes->numbermonitors = 0; 1326 PetscFunctionReturn(0); 1327 } 1328 1329 #undef __FUNCT__ 1330 #define __FUNCT__ "SNESSetConvergenceTest" 1331 /*@C 1332 SNESSetConvergenceTest - Sets the function that is to be used 1333 to test for convergence of the nonlinear iterative solution. 1334 1335 Collective on SNES 1336 1337 Input Parameters: 1338 + snes - the SNES context 1339 . func - routine to test for convergence 1340 - cctx - [optional] context for private data for the convergence routine 1341 (may be PETSC_NULL) 1342 1343 Calling sequence of func: 1344 $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1345 1346 + snes - the SNES context 1347 . cctx - [optional] convergence context 1348 . reason - reason for convergence/divergence 1349 . xnorm - 2-norm of current iterate 1350 . gnorm - 2-norm of current step 1351 - f - 2-norm of function 1352 1353 Level: advanced 1354 1355 .keywords: SNES, nonlinear, set, convergence, test 1356 1357 .seealso: SNESConverged_LS(), SNESConverged_TR() 1358 @*/ 1359 PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 1360 { 1361 PetscFunctionBegin; 1362 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1363 (snes)->converged = func; 1364 (snes)->cnvP = cctx; 1365 PetscFunctionReturn(0); 1366 } 1367 1368 #undef __FUNCT__ 1369 #define __FUNCT__ "SNESGetConvergedReason" 1370 /*@C 1371 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1372 1373 Not Collective 1374 1375 Input Parameter: 1376 . snes - the SNES context 1377 1378 Output Parameter: 1379 . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1380 manual pages for the individual convergence tests for complete lists 1381 1382 Level: intermediate 1383 1384 Notes: Can only be called after the call the SNESSolve() is complete. 1385 1386 .keywords: SNES, nonlinear, set, convergence, test 1387 1388 .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1389 @*/ 1390 PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1391 { 1392 PetscFunctionBegin; 1393 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1394 PetscValidPointer(reason,2); 1395 *reason = snes->reason; 1396 PetscFunctionReturn(0); 1397 } 1398 1399 #undef __FUNCT__ 1400 #define __FUNCT__ "SNESSetConvergenceHistory" 1401 /*@ 1402 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1403 1404 Collective on SNES 1405 1406 Input Parameters: 1407 + snes - iterative context obtained from SNESCreate() 1408 . a - array to hold history 1409 . its - integer array holds the number of linear iterations for each solve. 1410 . na - size of a and its 1411 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1412 else it continues storing new values for new nonlinear solves after the old ones 1413 1414 Notes: 1415 If set, this array will contain the function norms computed 1416 at each step. 1417 1418 This routine is useful, e.g., when running a code for purposes 1419 of accurate performance monitoring, when no I/O should be done 1420 during the section of code that is being timed. 1421 1422 Level: intermediate 1423 1424 .keywords: SNES, set, convergence, history 1425 1426 .seealso: SNESGetConvergenceHistory() 1427 1428 @*/ 1429 PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1430 { 1431 PetscFunctionBegin; 1432 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1433 if (na) PetscValidScalarPointer(a,2); 1434 snes->conv_hist = a; 1435 snes->conv_hist_its = its; 1436 snes->conv_hist_max = na; 1437 snes->conv_hist_reset = reset; 1438 PetscFunctionReturn(0); 1439 } 1440 1441 #undef __FUNCT__ 1442 #define __FUNCT__ "SNESGetConvergenceHistory" 1443 /*@C 1444 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1445 1446 Collective on SNES 1447 1448 Input Parameter: 1449 . snes - iterative context obtained from SNESCreate() 1450 1451 Output Parameters: 1452 . a - array to hold history 1453 . its - integer array holds the number of linear iterations (or 1454 negative if not converged) for each solve. 1455 - na - size of a and its 1456 1457 Notes: 1458 The calling sequence for this routine in Fortran is 1459 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1460 1461 This routine is useful, e.g., when running a code for purposes 1462 of accurate performance monitoring, when no I/O should be done 1463 during the section of code that is being timed. 1464 1465 Level: intermediate 1466 1467 .keywords: SNES, get, convergence, history 1468 1469 .seealso: SNESSetConvergencHistory() 1470 1471 @*/ 1472 PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1473 { 1474 PetscFunctionBegin; 1475 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1476 if (a) *a = snes->conv_hist; 1477 if (its) *its = snes->conv_hist_its; 1478 if (na) *na = snes->conv_hist_len; 1479 PetscFunctionReturn(0); 1480 } 1481 1482 #undef __FUNCT__ 1483 #define __FUNCT__ "SNESSetRhsBC" 1484 /*@C 1485 SNESSetRhsBC - Sets the function which applies boundary conditions 1486 to the Rhs of each system. 1487 1488 Collective on SNES 1489 1490 Input Parameters: 1491 . snes - The nonlinear solver context 1492 . func - The function 1493 1494 Calling sequence of func: 1495 . func (SNES snes, Vec rhs, void *ctx); 1496 1497 . rhs - The current rhs vector 1498 . ctx - The user-context 1499 1500 Level: intermediate 1501 1502 .keywords: SNES, Rhs, boundary conditions 1503 .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate() 1504 @*/ 1505 PetscErrorCode SNESSetRhsBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 1506 { 1507 PetscFunctionBegin; 1508 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1509 snes->applyrhsbc = func; 1510 PetscFunctionReturn(0); 1511 } 1512 1513 #undef __FUNCT__ 1514 #define __FUNCT__ "SNESDefaultRhsBC" 1515 /*@ 1516 SNESDefaultRhsBC - The default boundary condition function which does nothing. 1517 1518 Not collective 1519 1520 Input Parameters: 1521 . snes - The nonlinear solver context 1522 . rhs - The Rhs 1523 . ctx - The user-context 1524 1525 Level: beginner 1526 1527 .keywords: SNES, Rhs, boundary conditions 1528 .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate() 1529 @*/ 1530 PetscErrorCode SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx) 1531 { 1532 PetscFunctionBegin; 1533 PetscFunctionReturn(0); 1534 } 1535 1536 #undef __FUNCT__ 1537 #define __FUNCT__ "SNESSetSolutionBC" 1538 /*@C 1539 SNESSetSolutionBC - Sets the function which applies boundary conditions 1540 to the solution of each system. 1541 1542 Collective on SNES 1543 1544 Input Parameters: 1545 . snes - The nonlinear solver context 1546 . func - The function 1547 1548 Calling sequence of func: 1549 . func (SNES snes, Vec rsol, void *ctx); 1550 1551 . sol - The current solution vector 1552 . ctx - The user-context 1553 1554 Level: intermediate 1555 1556 .keywords: SNES, solution, boundary conditions 1557 .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate() 1558 @*/ 1559 PetscErrorCode SNESSetSolutionBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 1560 { 1561 PetscFunctionBegin; 1562 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1563 snes->applysolbc = func; 1564 PetscFunctionReturn(0); 1565 } 1566 1567 #undef __FUNCT__ 1568 #define __FUNCT__ "SNESDefaultSolutionBC" 1569 /*@ 1570 SNESDefaultSolutionBC - The default boundary condition function which does nothing. 1571 1572 Not collective 1573 1574 Input Parameters: 1575 . snes - The nonlinear solver context 1576 . sol - The solution 1577 . ctx - The user-context 1578 1579 Level: beginner 1580 1581 .keywords: SNES, solution, boundary conditions 1582 .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate() 1583 @*/ 1584 PetscErrorCode SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx) 1585 { 1586 PetscFunctionBegin; 1587 PetscFunctionReturn(0); 1588 } 1589 1590 #undef __FUNCT__ 1591 #define __FUNCT__ "SNESSetUpdate" 1592 /*@C 1593 SNESSetUpdate - Sets the general-purpose update function called 1594 at the beginning of every step of the iteration. 1595 1596 Collective on SNES 1597 1598 Input Parameters: 1599 . snes - The nonlinear solver context 1600 . func - The function 1601 1602 Calling sequence of func: 1603 . func (TS ts, int step); 1604 1605 . step - The current step of the iteration 1606 1607 Level: intermediate 1608 1609 .keywords: SNES, update 1610 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 1611 @*/ 1612 PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 1613 { 1614 PetscFunctionBegin; 1615 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1616 snes->update = func; 1617 PetscFunctionReturn(0); 1618 } 1619 1620 #undef __FUNCT__ 1621 #define __FUNCT__ "SNESDefaultUpdate" 1622 /*@ 1623 SNESDefaultUpdate - The default update function which does nothing. 1624 1625 Not collective 1626 1627 Input Parameters: 1628 . snes - The nonlinear solver context 1629 . step - The current step of the iteration 1630 1631 Level: intermediate 1632 1633 .keywords: SNES, update 1634 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 1635 @*/ 1636 PetscErrorCode SNESDefaultUpdate(SNES snes, PetscInt step) 1637 { 1638 PetscFunctionBegin; 1639 PetscFunctionReturn(0); 1640 } 1641 1642 #undef __FUNCT__ 1643 #define __FUNCT__ "SNESScaleStep_Private" 1644 /* 1645 SNESScaleStep_Private - Scales a step so that its length is less than the 1646 positive parameter delta. 1647 1648 Input Parameters: 1649 + snes - the SNES context 1650 . y - approximate solution of linear system 1651 . fnorm - 2-norm of current function 1652 - delta - trust region size 1653 1654 Output Parameters: 1655 + gpnorm - predicted function norm at the new point, assuming local 1656 linearization. The value is zero if the step lies within the trust 1657 region, and exceeds zero otherwise. 1658 - ynorm - 2-norm of the step 1659 1660 Note: 1661 For non-trust region methods such as SNESLS, the parameter delta 1662 is set to be the maximum allowable step size. 1663 1664 .keywords: SNES, nonlinear, scale, step 1665 */ 1666 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 1667 { 1668 PetscReal nrm; 1669 PetscScalar cnorm; 1670 PetscErrorCode ierr; 1671 1672 PetscFunctionBegin; 1673 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1674 PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1675 PetscCheckSameComm(snes,1,y,2); 1676 1677 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1678 if (nrm > *delta) { 1679 nrm = *delta/nrm; 1680 *gpnorm = (1.0 - nrm)*(*fnorm); 1681 cnorm = nrm; 1682 ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 1683 *ynorm = *delta; 1684 } else { 1685 *gpnorm = 0.0; 1686 *ynorm = nrm; 1687 } 1688 PetscFunctionReturn(0); 1689 } 1690 1691 static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero", 1692 "not currently used", 1693 "line search failed", 1694 "reach maximum number of iterations", 1695 "function norm became NaN (not a number)", 1696 "not currently used", 1697 "number of function computations exceeded", 1698 "not currently used", 1699 "still iterating", 1700 "not currently used", 1701 "absolute size of function norm", 1702 "relative decrease in function norm", 1703 "step size is small", 1704 "not currently used", 1705 "not currently used", 1706 "small trust region"}; 1707 1708 #undef __FUNCT__ 1709 #define __FUNCT__ "SNESSolve" 1710 /*@ 1711 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 1712 SNESCreate() and optional routines of the form SNESSetXXX(). 1713 1714 Collective on SNES 1715 1716 Input Parameters: 1717 + snes - the SNES context 1718 - x - the solution vector 1719 1720 Notes: 1721 The user should initialize the vector,x, with the initial guess 1722 for the nonlinear solve prior to calling SNESSolve. In particular, 1723 to employ an initial guess of zero, the user should explicitly set 1724 this vector to zero by calling VecSet(). 1725 1726 Level: beginner 1727 1728 .keywords: SNES, nonlinear, solve 1729 1730 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs() 1731 @*/ 1732 PetscErrorCode SNESSolve(SNES snes,Vec x) 1733 { 1734 PetscErrorCode ierr; 1735 PetscTruth flg; 1736 1737 PetscFunctionBegin; 1738 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1739 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1740 PetscCheckSameComm(snes,1,x,2); 1741 if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 1742 1743 if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1744 else {snes->vec_sol = snes->vec_sol_always = x;} 1745 if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 1746 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1747 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1748 ierr = (*(snes)->solve)(snes);CHKERRQ(ierr); 1749 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1750 ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 1751 if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1752 ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1753 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 1754 if (snes->printreason) { 1755 if (snes->reason > 0) { 1756 ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 1757 } else { 1758 ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 1759 } 1760 } 1761 1762 PetscFunctionReturn(0); 1763 } 1764 1765 /* --------- Internal routines for SNES Package --------- */ 1766 1767 #undef __FUNCT__ 1768 #define __FUNCT__ "SNESSetType" 1769 /*@C 1770 SNESSetType - Sets the method for the nonlinear solver. 1771 1772 Collective on SNES 1773 1774 Input Parameters: 1775 + snes - the SNES context 1776 - type - a known method 1777 1778 Options Database Key: 1779 . -snes_type <type> - Sets the method; use -help for a list 1780 of available methods (for instance, ls or tr) 1781 1782 Notes: 1783 See "petsc/include/petscsnes.h" for available methods (for instance) 1784 + SNESLS - Newton's method with line search 1785 (systems of nonlinear equations) 1786 . SNESTR - Newton's method with trust region 1787 (systems of nonlinear equations) 1788 1789 Normally, it is best to use the SNESSetFromOptions() command and then 1790 set the SNES solver type from the options database rather than by using 1791 this routine. Using the options database provides the user with 1792 maximum flexibility in evaluating the many nonlinear solvers. 1793 The SNESSetType() routine is provided for those situations where it 1794 is necessary to set the nonlinear solver independently of the command 1795 line or options database. This might be the case, for example, when 1796 the choice of solver changes during the execution of the program, 1797 and the user's application is taking responsibility for choosing the 1798 appropriate method. 1799 1800 Level: intermediate 1801 1802 .keywords: SNES, set, type 1803 1804 .seealso: SNESType, SNESCreate() 1805 1806 @*/ 1807 PetscErrorCode SNESSetType(SNES snes,const SNESType type) 1808 { 1809 PetscErrorCode ierr,(*r)(SNES); 1810 PetscTruth match; 1811 1812 PetscFunctionBegin; 1813 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1814 PetscValidCharPointer(type,2); 1815 1816 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 1817 if (match) PetscFunctionReturn(0); 1818 1819 if (snes->setupcalled) { 1820 ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 1821 snes->data = 0; 1822 } 1823 1824 /* Get the function pointers for the iterative method requested */ 1825 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1826 ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1827 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1828 if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 1829 snes->data = 0; 1830 ierr = (*r)(snes);CHKERRQ(ierr); 1831 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 1832 PetscFunctionReturn(0); 1833 } 1834 1835 1836 /* --------------------------------------------------------------------- */ 1837 #undef __FUNCT__ 1838 #define __FUNCT__ "SNESRegisterDestroy" 1839 /*@C 1840 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1841 registered by SNESRegisterDynamic(). 1842 1843 Not Collective 1844 1845 Level: advanced 1846 1847 .keywords: SNES, nonlinear, register, destroy 1848 1849 .seealso: SNESRegisterAll(), SNESRegisterAll() 1850 @*/ 1851 PetscErrorCode SNESRegisterDestroy(void) 1852 { 1853 PetscErrorCode ierr; 1854 1855 PetscFunctionBegin; 1856 if (SNESList) { 1857 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 1858 SNESList = 0; 1859 } 1860 SNESRegisterAllCalled = PETSC_FALSE; 1861 PetscFunctionReturn(0); 1862 } 1863 1864 #undef __FUNCT__ 1865 #define __FUNCT__ "SNESGetType" 1866 /*@C 1867 SNESGetType - Gets the SNES method type and name (as a string). 1868 1869 Not Collective 1870 1871 Input Parameter: 1872 . snes - nonlinear solver context 1873 1874 Output Parameter: 1875 . type - SNES method (a character string) 1876 1877 Level: intermediate 1878 1879 .keywords: SNES, nonlinear, get, type, name 1880 @*/ 1881 PetscErrorCode SNESGetType(SNES snes,SNESType *type) 1882 { 1883 PetscFunctionBegin; 1884 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1885 PetscValidPointer(type,2); 1886 *type = snes->type_name; 1887 PetscFunctionReturn(0); 1888 } 1889 1890 #undef __FUNCT__ 1891 #define __FUNCT__ "SNESGetSolution" 1892 /*@C 1893 SNESGetSolution - Returns the vector where the approximate solution is 1894 stored. 1895 1896 Not Collective, but Vec is parallel if SNES is parallel 1897 1898 Input Parameter: 1899 . snes - the SNES context 1900 1901 Output Parameter: 1902 . x - the solution 1903 1904 Level: advanced 1905 1906 .keywords: SNES, nonlinear, get, solution 1907 1908 .seealso: SNESGetFunction(), SNESGetSolutionUpdate() 1909 @*/ 1910 PetscErrorCode SNESGetSolution(SNES snes,Vec *x) 1911 { 1912 PetscFunctionBegin; 1913 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1914 PetscValidPointer(x,2); 1915 *x = snes->vec_sol_always; 1916 PetscFunctionReturn(0); 1917 } 1918 1919 #undef __FUNCT__ 1920 #define __FUNCT__ "SNESGetSolutionUpdate" 1921 /*@C 1922 SNESGetSolutionUpdate - Returns the vector where the solution update is 1923 stored. 1924 1925 Not Collective, but Vec is parallel if SNES is parallel 1926 1927 Input Parameter: 1928 . snes - the SNES context 1929 1930 Output Parameter: 1931 . x - the solution update 1932 1933 Level: advanced 1934 1935 .keywords: SNES, nonlinear, get, solution, update 1936 1937 .seealso: SNESGetSolution(), SNESGetFunction 1938 @*/ 1939 PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x) 1940 { 1941 PetscFunctionBegin; 1942 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1943 PetscValidPointer(x,2); 1944 *x = snes->vec_sol_update_always; 1945 PetscFunctionReturn(0); 1946 } 1947 1948 #undef __FUNCT__ 1949 #define __FUNCT__ "SNESGetFunction" 1950 /*@C 1951 SNESGetFunction - Returns the vector where the function is stored. 1952 1953 Not Collective, but Vec is parallel if SNES is parallel 1954 1955 Input Parameter: 1956 . snes - the SNES context 1957 1958 Output Parameter: 1959 + r - the function (or PETSC_NULL) 1960 . ctx - the function context (or PETSC_NULL) 1961 - func - the function (or PETSC_NULL) 1962 1963 Level: advanced 1964 1965 .keywords: SNES, nonlinear, get, function 1966 1967 .seealso: SNESSetFunction(), SNESGetSolution() 1968 @*/ 1969 PetscErrorCode SNESGetFunction(SNES snes,Vec *r,void **ctx,PetscErrorCode (**func)(SNES,Vec,Vec,void*)) 1970 { 1971 PetscFunctionBegin; 1972 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1973 if (r) *r = snes->vec_func_always; 1974 if (ctx) *ctx = snes->funP; 1975 if (func) *func = snes->computefunction; 1976 PetscFunctionReturn(0); 1977 } 1978 1979 #undef __FUNCT__ 1980 #define __FUNCT__ "SNESSetOptionsPrefix" 1981 /*@C 1982 SNESSetOptionsPrefix - Sets the prefix used for searching for all 1983 SNES options in the database. 1984 1985 Collective on SNES 1986 1987 Input Parameter: 1988 + snes - the SNES context 1989 - prefix - the prefix to prepend to all option names 1990 1991 Notes: 1992 A hyphen (-) must NOT be given at the beginning of the prefix name. 1993 The first character of all runtime options is AUTOMATICALLY the hyphen. 1994 1995 Level: advanced 1996 1997 .keywords: SNES, set, options, prefix, database 1998 1999 .seealso: SNESSetFromOptions() 2000 @*/ 2001 PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2002 { 2003 PetscErrorCode ierr; 2004 2005 PetscFunctionBegin; 2006 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2007 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2008 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2009 PetscFunctionReturn(0); 2010 } 2011 2012 #undef __FUNCT__ 2013 #define __FUNCT__ "SNESAppendOptionsPrefix" 2014 /*@C 2015 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2016 SNES options in the database. 2017 2018 Collective on SNES 2019 2020 Input Parameters: 2021 + snes - the SNES context 2022 - prefix - the prefix to prepend to all option names 2023 2024 Notes: 2025 A hyphen (-) must NOT be given at the beginning of the prefix name. 2026 The first character of all runtime options is AUTOMATICALLY the hyphen. 2027 2028 Level: advanced 2029 2030 .keywords: SNES, append, options, prefix, database 2031 2032 .seealso: SNESGetOptionsPrefix() 2033 @*/ 2034 PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2035 { 2036 PetscErrorCode ierr; 2037 2038 PetscFunctionBegin; 2039 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2040 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2041 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2042 PetscFunctionReturn(0); 2043 } 2044 2045 #undef __FUNCT__ 2046 #define __FUNCT__ "SNESGetOptionsPrefix" 2047 /*@C 2048 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2049 SNES options in the database. 2050 2051 Not Collective 2052 2053 Input Parameter: 2054 . snes - the SNES context 2055 2056 Output Parameter: 2057 . prefix - pointer to the prefix string used 2058 2059 Notes: On the fortran side, the user should pass in a string 'prifix' of 2060 sufficient length to hold the prefix. 2061 2062 Level: advanced 2063 2064 .keywords: SNES, get, options, prefix, database 2065 2066 .seealso: SNESAppendOptionsPrefix() 2067 @*/ 2068 PetscErrorCode SNESGetOptionsPrefix(SNES snes,char *prefix[]) 2069 { 2070 PetscErrorCode ierr; 2071 2072 PetscFunctionBegin; 2073 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2074 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2075 PetscFunctionReturn(0); 2076 } 2077 2078 2079 #undef __FUNCT__ 2080 #define __FUNCT__ "SNESRegister" 2081 /*@C 2082 SNESRegister - See SNESRegisterDynamic() 2083 2084 Level: advanced 2085 @*/ 2086 PetscErrorCode SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2087 { 2088 char fullname[PETSC_MAX_PATH_LEN]; 2089 PetscErrorCode ierr; 2090 2091 PetscFunctionBegin; 2092 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2093 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2094 PetscFunctionReturn(0); 2095 } 2096 2097 #undef __FUNCT__ 2098 #define __FUNCT__ "SNESTestLocalMin" 2099 PetscErrorCode SNESTestLocalMin(SNES snes) 2100 { 2101 PetscErrorCode ierr; 2102 PetscInt N,i,j; 2103 Vec u,uh,fh; 2104 PetscScalar value; 2105 PetscReal norm; 2106 2107 PetscFunctionBegin; 2108 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2109 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2110 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2111 2112 /* currently only works for sequential */ 2113 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2114 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2115 for (i=0; i<N; i++) { 2116 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2117 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2118 for (j=-10; j<11; j++) { 2119 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2120 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2121 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2122 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2123 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2124 value = -value; 2125 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2126 } 2127 } 2128 ierr = VecDestroy(uh);CHKERRQ(ierr); 2129 ierr = VecDestroy(fh);CHKERRQ(ierr); 2130 PetscFunctionReturn(0); 2131 } 2132