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 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 snes->bops->publish = SNESPublish_Petsc; 595 snes->max_its = 50; 596 snes->max_funcs = 10000; 597 snes->norm = 0.0; 598 snes->rtol = 1.e-8; 599 snes->ttol = 0.0; 600 snes->abstol = 1.e-50; 601 snes->xtol = 1.e-8; 602 snes->deltatol = 1.e-12; 603 snes->nfuncs = 0; 604 snes->numFailures = 0; 605 snes->maxFailures = 1; 606 snes->linear_its = 0; 607 snes->numbermonitors = 0; 608 snes->data = 0; 609 snes->view = 0; 610 snes->setupcalled = 0; 611 snes->ksp_ewconv = PETSC_FALSE; 612 snes->vwork = 0; 613 snes->nwork = 0; 614 snes->conv_hist_len = 0; 615 snes->conv_hist_max = 0; 616 snes->conv_hist = PETSC_NULL; 617 snes->conv_hist_its = PETSC_NULL; 618 snes->conv_hist_reset = PETSC_TRUE; 619 snes->reason = SNES_CONVERGED_ITERATING; 620 621 /* Create context to compute Eisenstat-Walker relative tolerance for KSP */ 622 ierr = PetscNew(SNES_KSP_EW_ConvCtx,&kctx);CHKERRQ(ierr); 623 PetscLogObjectMemory(snes,sizeof(SNES_KSP_EW_ConvCtx)); 624 snes->kspconvctx = (void*)kctx; 625 kctx->version = 2; 626 kctx->rtol_0 = .3; /* Eisenstat and Walker suggest rtol_0=.5, but 627 this was too large for some test cases */ 628 kctx->rtol_last = 0; 629 kctx->rtol_max = .9; 630 kctx->gamma = 1.0; 631 kctx->alpha2 = .5*(1.0 + sqrt(5.0)); 632 kctx->alpha = kctx->alpha2; 633 kctx->threshold = .1; 634 kctx->lresid_last = 0; 635 kctx->norm_last = 0; 636 637 ierr = KSPCreate(comm,&snes->ksp);CHKERRQ(ierr); 638 PetscLogObjectParent(snes,snes->ksp) 639 640 *outsnes = snes; 641 ierr = PetscPublishAll(snes);CHKERRQ(ierr); 642 PetscFunctionReturn(0); 643 } 644 645 #undef __FUNCT__ 646 #define __FUNCT__ "SNESSetFunction" 647 /*@C 648 SNESSetFunction - Sets the function evaluation routine and function 649 vector for use by the SNES routines in solving systems of nonlinear 650 equations. 651 652 Collective on SNES 653 654 Input Parameters: 655 + snes - the SNES context 656 . func - function evaluation routine 657 . r - vector to store function value 658 - ctx - [optional] user-defined context for private data for the 659 function evaluation routine (may be PETSC_NULL) 660 661 Calling sequence of func: 662 $ func (SNES snes,Vec x,Vec f,void *ctx); 663 664 . f - function vector 665 - ctx - optional user-defined function context 666 667 Notes: 668 The Newton-like methods typically solve linear systems of the form 669 $ f'(x) x = -f(x), 670 where f'(x) denotes the Jacobian matrix and f(x) is the function. 671 672 Level: beginner 673 674 .keywords: SNES, nonlinear, set, function 675 676 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian() 677 @*/ 678 PetscErrorCode SNESSetFunction(SNES snes,Vec r,PetscErrorCode (*func)(SNES,Vec,Vec,void*),void *ctx) 679 { 680 PetscFunctionBegin; 681 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 682 PetscValidHeaderSpecific(r,VEC_COOKIE,2); 683 PetscCheckSameComm(snes,1,r,2); 684 685 snes->computefunction = func; 686 snes->vec_func = snes->vec_func_always = r; 687 snes->funP = ctx; 688 PetscFunctionReturn(0); 689 } 690 691 /* --------------------------------------------------------------- */ 692 #undef __FUNCT__ 693 #define __FUNCT__ "SNESSetRhs" 694 /*@C 695 SNESSetRhs - Sets the vector for solving F(x) = rhs. If rhs is not set 696 it assumes a zero right hand side. 697 698 Collective on SNES 699 700 Input Parameters: 701 + snes - the SNES context 702 - rhs - the right hand side vector or PETSC_NULL for a zero right hand side 703 704 Level: intermediate 705 706 .keywords: SNES, nonlinear, set, function, right hand side 707 708 .seealso: SNESGetFunction(), SNESComputeFunction(), SNESSetJacobian(), SNESSetFunction() 709 @*/ 710 PetscErrorCode SNESSetRhs(SNES snes,Vec rhs) 711 { 712 PetscErrorCode ierr; 713 714 PetscFunctionBegin; 715 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 716 if (rhs) { 717 PetscValidHeaderSpecific(rhs,VEC_COOKIE,2); 718 PetscCheckSameComm(snes,1,rhs,2); 719 ierr = PetscObjectReference((PetscObject)rhs);CHKERRQ(ierr); 720 } 721 if (snes->afine) { 722 ierr = VecDestroy(snes->afine);CHKERRQ(ierr); 723 } 724 snes->afine = rhs; 725 PetscFunctionReturn(0); 726 } 727 728 #undef __FUNCT__ 729 #define __FUNCT__ "SNESComputeFunction" 730 /*@ 731 SNESComputeFunction - Calls the function that has been set with 732 SNESSetFunction(). 733 734 Collective on SNES 735 736 Input Parameters: 737 + snes - the SNES context 738 - x - input vector 739 740 Output Parameter: 741 . y - function vector, as set by SNESSetFunction() 742 743 Notes: 744 SNESComputeFunction() is typically used within nonlinear solvers 745 implementations, so most users would not generally call this routine 746 themselves. 747 748 Level: developer 749 750 .keywords: SNES, nonlinear, compute, function 751 752 .seealso: SNESSetFunction(), SNESGetFunction() 753 @*/ 754 PetscErrorCode SNESComputeFunction(SNES snes,Vec x,Vec y) 755 { 756 PetscErrorCode ierr; 757 758 PetscFunctionBegin; 759 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 760 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 761 PetscValidHeaderSpecific(y,VEC_COOKIE,3); 762 PetscCheckSameComm(snes,1,x,2); 763 PetscCheckSameComm(snes,1,y,3); 764 765 ierr = PetscLogEventBegin(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 766 PetscStackPush("SNES user function"); 767 ierr = (*snes->computefunction)(snes,x,y,snes->funP);CHKERRQ(ierr); 768 PetscStackPop; 769 if (snes->afine) { 770 PetscScalar mone = -1.0; 771 ierr = VecAXPY(&mone,snes->afine,y);CHKERRQ(ierr); 772 } 773 snes->nfuncs++; 774 ierr = PetscLogEventEnd(SNES_FunctionEval,snes,x,y,0);CHKERRQ(ierr); 775 PetscFunctionReturn(0); 776 } 777 778 #undef __FUNCT__ 779 #define __FUNCT__ "SNESComputeJacobian" 780 /*@ 781 SNESComputeJacobian - Computes the Jacobian matrix that has been 782 set with SNESSetJacobian(). 783 784 Collective on SNES and Mat 785 786 Input Parameters: 787 + snes - the SNES context 788 - x - input vector 789 790 Output Parameters: 791 + A - Jacobian matrix 792 . B - optional preconditioning matrix 793 - flag - flag indicating matrix structure 794 795 Notes: 796 Most users should not need to explicitly call this routine, as it 797 is used internally within the nonlinear solvers. 798 799 See KSPSetOperators() for important information about setting the 800 flag parameter. 801 802 Level: developer 803 804 .keywords: SNES, compute, Jacobian, matrix 805 806 .seealso: SNESSetJacobian(), KSPSetOperators() 807 @*/ 808 PetscErrorCode SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 809 { 810 PetscErrorCode ierr; 811 812 PetscFunctionBegin; 813 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 814 PetscValidHeaderSpecific(X,VEC_COOKIE,2); 815 PetscValidPointer(flg,5); 816 PetscCheckSameComm(snes,1,X,2); 817 if (!snes->computejacobian) PetscFunctionReturn(0); 818 ierr = PetscLogEventBegin(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 819 *flg = DIFFERENT_NONZERO_PATTERN; 820 PetscStackPush("SNES user Jacobian function"); 821 ierr = (*snes->computejacobian)(snes,X,A,B,flg,snes->jacP);CHKERRQ(ierr); 822 PetscStackPop; 823 ierr = PetscLogEventEnd(SNES_JacobianEval,snes,X,*A,*B);CHKERRQ(ierr); 824 /* make sure user returned a correct Jacobian and preconditioner */ 825 PetscValidHeaderSpecific(*A,MAT_COOKIE,3); 826 PetscValidHeaderSpecific(*B,MAT_COOKIE,4); 827 PetscFunctionReturn(0); 828 } 829 830 #undef __FUNCT__ 831 #define __FUNCT__ "SNESSetJacobian" 832 /*@C 833 SNESSetJacobian - Sets the function to compute Jacobian as well as the 834 location to store the matrix. 835 836 Collective on SNES and Mat 837 838 Input Parameters: 839 + snes - the SNES context 840 . A - Jacobian matrix 841 . B - preconditioner matrix (usually same as the Jacobian) 842 . func - Jacobian evaluation routine 843 - ctx - [optional] user-defined context for private data for the 844 Jacobian evaluation routine (may be PETSC_NULL) 845 846 Calling sequence of func: 847 $ func (SNES snes,Vec x,Mat *A,Mat *B,int *flag,void *ctx); 848 849 + x - input vector 850 . A - Jacobian matrix 851 . B - preconditioner matrix, usually the same as A 852 . flag - flag indicating information about the preconditioner matrix 853 structure (same as flag in KSPSetOperators()) 854 - ctx - [optional] user-defined Jacobian context 855 856 Notes: 857 See KSPSetOperators() for important information about setting the flag 858 output parameter in the routine func(). Be sure to read this information! 859 860 The routine func() takes Mat * as the matrix arguments rather than Mat. 861 This allows the Jacobian evaluation routine to replace A and/or B with a 862 completely new new matrix structure (not just different matrix elements) 863 when appropriate, for instance, if the nonzero structure is changing 864 throughout the global iterations. 865 866 Level: beginner 867 868 .keywords: SNES, nonlinear, set, Jacobian, matrix 869 870 .seealso: KSPSetOperators(), SNESSetFunction(), , MatSNESMFComputeJacobian(), SNESDefaultComputeJacobianColor() 871 @*/ 872 PetscErrorCode SNESSetJacobian(SNES snes,Mat A,Mat B,PetscErrorCode (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 873 { 874 PetscErrorCode ierr; 875 876 PetscFunctionBegin; 877 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 878 if (A) PetscValidHeaderSpecific(A,MAT_COOKIE,2); 879 if (B) PetscValidHeaderSpecific(B,MAT_COOKIE,3); 880 if (A) PetscCheckSameComm(snes,1,A,2); 881 if (B) PetscCheckSameComm(snes,1,B,2); 882 if (func) snes->computejacobian = func; 883 if (ctx) snes->jacP = ctx; 884 if (A) { 885 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 886 snes->jacobian = A; 887 ierr = PetscObjectReference((PetscObject)A);CHKERRQ(ierr); 888 } 889 if (B) { 890 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 891 snes->jacobian_pre = B; 892 ierr = PetscObjectReference((PetscObject)B);CHKERRQ(ierr); 893 } 894 ierr = KSPSetOperators(snes->ksp,A,B,SAME_NONZERO_PATTERN);CHKERRQ(ierr); 895 PetscFunctionReturn(0); 896 } 897 898 #undef __FUNCT__ 899 #define __FUNCT__ "SNESGetJacobian" 900 /*@C 901 SNESGetJacobian - Returns the Jacobian matrix and optionally the user 902 provided context for evaluating the Jacobian. 903 904 Not Collective, but Mat object will be parallel if SNES object is 905 906 Input Parameter: 907 . snes - the nonlinear solver context 908 909 Output Parameters: 910 + A - location to stash Jacobian matrix (or PETSC_NULL) 911 . B - location to stash preconditioner matrix (or PETSC_NULL) 912 . ctx - location to stash Jacobian ctx (or PETSC_NULL) 913 - func - location to put Jacobian function (or PETSC_NULL) 914 915 Level: advanced 916 917 .seealso: SNESSetJacobian(), SNESComputeJacobian() 918 @*/ 919 PetscErrorCode SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,PetscErrorCode (**func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*)) 920 { 921 PetscFunctionBegin; 922 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 923 if (A) *A = snes->jacobian; 924 if (B) *B = snes->jacobian_pre; 925 if (ctx) *ctx = snes->jacP; 926 if (func) *func = snes->computejacobian; 927 PetscFunctionReturn(0); 928 } 929 930 /* ----- Routines to initialize and destroy a nonlinear solver ---- */ 931 EXTERN PetscErrorCode SNESDefaultMatrixFreeCreate2(SNES,Vec,Mat*); 932 933 #undef __FUNCT__ 934 #define __FUNCT__ "SNESSetUp" 935 /*@ 936 SNESSetUp - Sets up the internal data structures for the later use 937 of a nonlinear solver. 938 939 Collective on SNES 940 941 Input Parameters: 942 + snes - the SNES context 943 - x - the solution vector 944 945 Notes: 946 For basic use of the SNES solvers the user need not explicitly call 947 SNESSetUp(), since these actions will automatically occur during 948 the call to SNESSolve(). However, if one wishes to control this 949 phase separately, SNESSetUp() should be called after SNESCreate() 950 and optional routines of the form SNESSetXXX(), but before SNESSolve(). 951 952 Level: advanced 953 954 .keywords: SNES, nonlinear, setup 955 956 .seealso: SNESCreate(), SNESSolve(), SNESDestroy() 957 @*/ 958 PetscErrorCode SNESSetUp(SNES snes,Vec x) 959 { 960 PetscErrorCode ierr; 961 PetscTruth flg, iseqtr; 962 963 PetscFunctionBegin; 964 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 965 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 966 PetscCheckSameComm(snes,1,x,2); 967 snes->vec_sol = snes->vec_sol_always = x; 968 969 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator",&flg);CHKERRQ(ierr); 970 /* 971 This version replaces the user provided Jacobian matrix with a 972 matrix-free version but still employs the user-provided preconditioner matrix 973 */ 974 if (flg) { 975 Mat J; 976 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 977 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 978 PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator routines\n"); 979 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 980 ierr = MatDestroy(J);CHKERRQ(ierr); 981 } 982 983 #if !defined(PETSC_USE_COMPLEX) && !defined(PETSC_USE_SINGLE) 984 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf_operator2",&flg);CHKERRQ(ierr); 985 if (flg) { 986 Mat J; 987 ierr = SNESDefaultMatrixFreeCreate2(snes,snes->vec_sol,&J);CHKERRQ(ierr); 988 ierr = SNESSetJacobian(snes,J,0,0,0);CHKERRQ(ierr); 989 ierr = MatDestroy(J);CHKERRQ(ierr); 990 } 991 #endif 992 993 ierr = PetscOptionsHasName(snes->prefix,"-snes_mf",&flg);CHKERRQ(ierr); 994 /* 995 This version replaces both the user-provided Jacobian and the user- 996 provided preconditioner matrix with the default matrix free version. 997 */ 998 if (flg) { 999 Mat J; 1000 KSP ksp; 1001 PC pc; 1002 1003 ierr = MatCreateSNESMF(snes,snes->vec_sol,&J);CHKERRQ(ierr); 1004 ierr = MatSNESMFSetFromOptions(J);CHKERRQ(ierr); 1005 PetscLogInfo(snes,"SNESSetUp: Setting default matrix-free operator and preconditioner routines\n"); 1006 ierr = SNESSetJacobian(snes,J,J,MatSNESMFComputeJacobian,snes->funP);CHKERRQ(ierr); 1007 ierr = MatDestroy(J);CHKERRQ(ierr); 1008 1009 /* force no preconditioner */ 1010 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1011 ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1012 ierr = PetscTypeCompare((PetscObject)pc,PCSHELL,&flg);CHKERRQ(ierr); 1013 if (!flg) { 1014 ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); 1015 } 1016 } 1017 1018 if (!snes->vec_func) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1019 if (!snes->computefunction) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetFunction() first"); 1020 if (!snes->jacobian) SETERRQ(PETSC_ERR_ARG_WRONGSTATE,"Must call SNESSetJacobian() first \n or use -snes_mf option"); 1021 if (snes->vec_func == snes->vec_sol) { 1022 SETERRQ(PETSC_ERR_ARG_IDN,"Solution vector cannot be function vector"); 1023 } 1024 1025 /* Set the KSP stopping criterion to use the Eisenstat-Walker method */ 1026 ierr = PetscTypeCompare((PetscObject)snes,SNESTR,&iseqtr);CHKERRQ(ierr); 1027 if (snes->ksp_ewconv && !iseqtr) { 1028 KSP ksp; 1029 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 1030 ierr = KSPSetConvergenceTest(ksp,SNES_KSP_EW_Converged_Private,snes);CHKERRQ(ierr); 1031 } 1032 1033 if (snes->setup) {ierr = (*snes->setup)(snes);CHKERRQ(ierr);} 1034 snes->setupcalled = 1; 1035 PetscFunctionReturn(0); 1036 } 1037 1038 #undef __FUNCT__ 1039 #define __FUNCT__ "SNESDestroy" 1040 /*@C 1041 SNESDestroy - Destroys the nonlinear solver context that was created 1042 with SNESCreate(). 1043 1044 Collective on SNES 1045 1046 Input Parameter: 1047 . snes - the SNES context 1048 1049 Level: beginner 1050 1051 .keywords: SNES, nonlinear, destroy 1052 1053 .seealso: SNESCreate(), SNESSolve() 1054 @*/ 1055 PetscErrorCode SNESDestroy(SNES snes) 1056 { 1057 PetscInt i; 1058 PetscErrorCode ierr; 1059 1060 PetscFunctionBegin; 1061 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1062 if (--snes->refct > 0) PetscFunctionReturn(0); 1063 1064 /* if memory was published with AMS then destroy it */ 1065 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1066 1067 if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1068 if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 1069 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1070 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1071 if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 1072 ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1073 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1074 for (i=0; i<snes->numbermonitors; i++) { 1075 if (snes->monitordestroy[i]) { 1076 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1077 } 1078 } 1079 ierr = PetscHeaderDestroy((PetscObject)snes);CHKERRQ(ierr); 1080 PetscFunctionReturn(0); 1081 } 1082 1083 /* ----------- Routines to set solver parameters ---------- */ 1084 1085 #undef __FUNCT__ 1086 #define __FUNCT__ "SNESSetTolerances" 1087 /*@ 1088 SNESSetTolerances - Sets various parameters used in convergence tests. 1089 1090 Collective on SNES 1091 1092 Input Parameters: 1093 + snes - the SNES context 1094 . abstol - absolute convergence tolerance 1095 . rtol - relative convergence tolerance 1096 . stol - convergence tolerance in terms of the norm 1097 of the change in the solution between steps 1098 . maxit - maximum number of iterations 1099 - maxf - maximum number of function evaluations 1100 1101 Options Database Keys: 1102 + -snes_atol <abstol> - Sets abstol 1103 . -snes_rtol <rtol> - Sets rtol 1104 . -snes_stol <stol> - Sets stol 1105 . -snes_max_it <maxit> - Sets maxit 1106 - -snes_max_funcs <maxf> - Sets maxf 1107 1108 Notes: 1109 The default maximum number of iterations is 50. 1110 The default maximum number of function evaluations is 1000. 1111 1112 Level: intermediate 1113 1114 .keywords: SNES, nonlinear, set, convergence, tolerances 1115 1116 .seealso: SNESSetTrustRegionTolerance() 1117 @*/ 1118 PetscErrorCode SNESSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit,PetscInt maxf) 1119 { 1120 PetscFunctionBegin; 1121 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1122 if (abstol != PETSC_DEFAULT) snes->abstol = abstol; 1123 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1124 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1125 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1126 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1127 PetscFunctionReturn(0); 1128 } 1129 1130 #undef __FUNCT__ 1131 #define __FUNCT__ "SNESGetTolerances" 1132 /*@ 1133 SNESGetTolerances - Gets various parameters used in convergence tests. 1134 1135 Not Collective 1136 1137 Input Parameters: 1138 + snes - the SNES context 1139 . abstol - absolute convergence tolerance 1140 . rtol - relative convergence tolerance 1141 . stol - convergence tolerance in terms of the norm 1142 of the change in the solution between steps 1143 . maxit - maximum number of iterations 1144 - maxf - maximum number of function evaluations 1145 1146 Notes: 1147 The user can specify PETSC_NULL for any parameter that is not needed. 1148 1149 Level: intermediate 1150 1151 .keywords: SNES, nonlinear, get, convergence, tolerances 1152 1153 .seealso: SNESSetTolerances() 1154 @*/ 1155 PetscErrorCode SNESGetTolerances(SNES snes,PetscReal *abstol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit,PetscInt *maxf) 1156 { 1157 PetscFunctionBegin; 1158 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1159 if (abstol) *abstol = snes->abstol; 1160 if (rtol) *rtol = snes->rtol; 1161 if (stol) *stol = snes->xtol; 1162 if (maxit) *maxit = snes->max_its; 1163 if (maxf) *maxf = snes->max_funcs; 1164 PetscFunctionReturn(0); 1165 } 1166 1167 #undef __FUNCT__ 1168 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1169 /*@ 1170 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1171 1172 Collective on SNES 1173 1174 Input Parameters: 1175 + snes - the SNES context 1176 - tol - tolerance 1177 1178 Options Database Key: 1179 . -snes_trtol <tol> - Sets tol 1180 1181 Level: intermediate 1182 1183 .keywords: SNES, nonlinear, set, trust region, tolerance 1184 1185 .seealso: SNESSetTolerances() 1186 @*/ 1187 PetscErrorCode SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1188 { 1189 PetscFunctionBegin; 1190 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1191 snes->deltatol = tol; 1192 PetscFunctionReturn(0); 1193 } 1194 1195 /* 1196 Duplicate the lg monitors for SNES from KSP; for some reason with 1197 dynamic libraries things don't work under Sun4 if we just use 1198 macros instead of functions 1199 */ 1200 #undef __FUNCT__ 1201 #define __FUNCT__ "SNESLGMonitor" 1202 PetscErrorCode SNESLGMonitor(SNES snes,PetscInt it,PetscReal norm,void *ctx) 1203 { 1204 PetscErrorCode ierr; 1205 1206 PetscFunctionBegin; 1207 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1208 ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1209 PetscFunctionReturn(0); 1210 } 1211 1212 #undef __FUNCT__ 1213 #define __FUNCT__ "SNESLGMonitorCreate" 1214 PetscErrorCode SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1215 { 1216 PetscErrorCode ierr; 1217 1218 PetscFunctionBegin; 1219 ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1220 PetscFunctionReturn(0); 1221 } 1222 1223 #undef __FUNCT__ 1224 #define __FUNCT__ "SNESLGMonitorDestroy" 1225 PetscErrorCode SNESLGMonitorDestroy(PetscDrawLG draw) 1226 { 1227 PetscErrorCode ierr; 1228 1229 PetscFunctionBegin; 1230 ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1231 PetscFunctionReturn(0); 1232 } 1233 1234 /* ------------ Routines to set performance monitoring options ----------- */ 1235 1236 #undef __FUNCT__ 1237 #define __FUNCT__ "SNESSetMonitor" 1238 /*@C 1239 SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 1240 iteration of the nonlinear solver to display the iteration's 1241 progress. 1242 1243 Collective on SNES 1244 1245 Input Parameters: 1246 + snes - the SNES context 1247 . func - monitoring routine 1248 . mctx - [optional] user-defined context for private data for the 1249 monitor routine (use PETSC_NULL if no context is desitre) 1250 - monitordestroy - [optional] routine that frees monitor context 1251 (may be PETSC_NULL) 1252 1253 Calling sequence of func: 1254 $ int func(SNES snes,PetscInt its, PetscReal norm,void *mctx) 1255 1256 + snes - the SNES context 1257 . its - iteration number 1258 . norm - 2-norm function value (may be estimated) 1259 - mctx - [optional] monitoring context 1260 1261 Options Database Keys: 1262 + -snes_monitor - sets SNESDefaultMonitor() 1263 . -snes_xmonitor - sets line graph monitor, 1264 uses SNESLGMonitorCreate() 1265 _ -snes_cancelmonitors - cancels all monitors that have 1266 been hardwired into a code by 1267 calls to SNESSetMonitor(), but 1268 does not cancel those set via 1269 the options database. 1270 1271 Notes: 1272 Several different monitoring routines may be set by calling 1273 SNESSetMonitor() multiple times; all will be called in the 1274 order in which they were set. 1275 1276 Level: intermediate 1277 1278 .keywords: SNES, nonlinear, set, monitor 1279 1280 .seealso: SNESDefaultMonitor(), SNESClearMonitor() 1281 @*/ 1282 PetscErrorCode SNESSetMonitor(SNES snes,PetscErrorCode (*func)(SNES,PetscInt,PetscReal,void*),void *mctx,PetscErrorCode (*monitordestroy)(void*)) 1283 { 1284 PetscFunctionBegin; 1285 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1286 if (snes->numbermonitors >= MAXSNESMONITORS) { 1287 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1288 } 1289 snes->monitor[snes->numbermonitors] = func; 1290 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1291 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1292 PetscFunctionReturn(0); 1293 } 1294 1295 #undef __FUNCT__ 1296 #define __FUNCT__ "SNESClearMonitor" 1297 /*@C 1298 SNESClearMonitor - Clears all the monitor functions for a SNES object. 1299 1300 Collective on SNES 1301 1302 Input Parameters: 1303 . snes - the SNES context 1304 1305 Options Database Key: 1306 . -snes_cancelmonitors - cancels all monitors that have been hardwired 1307 into a code by calls to SNESSetMonitor(), but does not cancel those 1308 set via the options database 1309 1310 Notes: 1311 There is no way to clear one specific monitor from a SNES object. 1312 1313 Level: intermediate 1314 1315 .keywords: SNES, nonlinear, set, monitor 1316 1317 .seealso: SNESDefaultMonitor(), SNESSetMonitor() 1318 @*/ 1319 PetscErrorCode SNESClearMonitor(SNES snes) 1320 { 1321 PetscFunctionBegin; 1322 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1323 snes->numbermonitors = 0; 1324 PetscFunctionReturn(0); 1325 } 1326 1327 #undef __FUNCT__ 1328 #define __FUNCT__ "SNESSetConvergenceTest" 1329 /*@C 1330 SNESSetConvergenceTest - Sets the function that is to be used 1331 to test for convergence of the nonlinear iterative solution. 1332 1333 Collective on SNES 1334 1335 Input Parameters: 1336 + snes - the SNES context 1337 . func - routine to test for convergence 1338 - cctx - [optional] context for private data for the convergence routine 1339 (may be PETSC_NULL) 1340 1341 Calling sequence of func: 1342 $ PetscErrorCode func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1343 1344 + snes - the SNES context 1345 . cctx - [optional] convergence context 1346 . reason - reason for convergence/divergence 1347 . xnorm - 2-norm of current iterate 1348 . gnorm - 2-norm of current step 1349 - f - 2-norm of function 1350 1351 Level: advanced 1352 1353 .keywords: SNES, nonlinear, set, convergence, test 1354 1355 .seealso: SNESConverged_LS(), SNESConverged_TR() 1356 @*/ 1357 PetscErrorCode SNESSetConvergenceTest(SNES snes,PetscErrorCode (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 1358 { 1359 PetscFunctionBegin; 1360 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1361 (snes)->converged = func; 1362 (snes)->cnvP = cctx; 1363 PetscFunctionReturn(0); 1364 } 1365 1366 #undef __FUNCT__ 1367 #define __FUNCT__ "SNESGetConvergedReason" 1368 /*@C 1369 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1370 1371 Not Collective 1372 1373 Input Parameter: 1374 . snes - the SNES context 1375 1376 Output Parameter: 1377 . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1378 manual pages for the individual convergence tests for complete lists 1379 1380 Level: intermediate 1381 1382 Notes: Can only be called after the call the SNESSolve() is complete. 1383 1384 .keywords: SNES, nonlinear, set, convergence, test 1385 1386 .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1387 @*/ 1388 PetscErrorCode SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1389 { 1390 PetscFunctionBegin; 1391 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1392 PetscValidPointer(reason,2); 1393 *reason = snes->reason; 1394 PetscFunctionReturn(0); 1395 } 1396 1397 #undef __FUNCT__ 1398 #define __FUNCT__ "SNESSetConvergenceHistory" 1399 /*@ 1400 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1401 1402 Collective on SNES 1403 1404 Input Parameters: 1405 + snes - iterative context obtained from SNESCreate() 1406 . a - array to hold history 1407 . its - integer array holds the number of linear iterations for each solve. 1408 . na - size of a and its 1409 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1410 else it continues storing new values for new nonlinear solves after the old ones 1411 1412 Notes: 1413 If set, this array will contain the function norms computed 1414 at each step. 1415 1416 This routine is useful, e.g., when running a code for purposes 1417 of accurate performance monitoring, when no I/O should be done 1418 during the section of code that is being timed. 1419 1420 Level: intermediate 1421 1422 .keywords: SNES, set, convergence, history 1423 1424 .seealso: SNESGetConvergenceHistory() 1425 1426 @*/ 1427 PetscErrorCode SNESSetConvergenceHistory(SNES snes,PetscReal a[],PetscInt *its,PetscInt na,PetscTruth reset) 1428 { 1429 PetscFunctionBegin; 1430 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1431 if (na) PetscValidScalarPointer(a,2); 1432 snes->conv_hist = a; 1433 snes->conv_hist_its = its; 1434 snes->conv_hist_max = na; 1435 snes->conv_hist_reset = reset; 1436 PetscFunctionReturn(0); 1437 } 1438 1439 #undef __FUNCT__ 1440 #define __FUNCT__ "SNESGetConvergenceHistory" 1441 /*@C 1442 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1443 1444 Collective on SNES 1445 1446 Input Parameter: 1447 . snes - iterative context obtained from SNESCreate() 1448 1449 Output Parameters: 1450 . a - array to hold history 1451 . its - integer array holds the number of linear iterations (or 1452 negative if not converged) for each solve. 1453 - na - size of a and its 1454 1455 Notes: 1456 The calling sequence for this routine in Fortran is 1457 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1458 1459 This routine is useful, e.g., when running a code for purposes 1460 of accurate performance monitoring, when no I/O should be done 1461 during the section of code that is being timed. 1462 1463 Level: intermediate 1464 1465 .keywords: SNES, get, convergence, history 1466 1467 .seealso: SNESSetConvergencHistory() 1468 1469 @*/ 1470 PetscErrorCode SNESGetConvergenceHistory(SNES snes,PetscReal *a[],PetscInt *its[],PetscInt *na) 1471 { 1472 PetscFunctionBegin; 1473 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1474 if (a) *a = snes->conv_hist; 1475 if (its) *its = snes->conv_hist_its; 1476 if (na) *na = snes->conv_hist_len; 1477 PetscFunctionReturn(0); 1478 } 1479 1480 #undef __FUNCT__ 1481 #define __FUNCT__ "SNESSetRhsBC" 1482 /*@C 1483 SNESSetRhsBC - Sets the function which applies boundary conditions 1484 to the Rhs of each system. 1485 1486 Collective on SNES 1487 1488 Input Parameters: 1489 . snes - The nonlinear solver context 1490 . func - The function 1491 1492 Calling sequence of func: 1493 . func (SNES snes, Vec rhs, void *ctx); 1494 1495 . rhs - The current rhs vector 1496 . ctx - The user-context 1497 1498 Level: intermediate 1499 1500 .keywords: SNES, Rhs, boundary conditions 1501 .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate() 1502 @*/ 1503 PetscErrorCode SNESSetRhsBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 1504 { 1505 PetscFunctionBegin; 1506 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1507 snes->applyrhsbc = func; 1508 PetscFunctionReturn(0); 1509 } 1510 1511 #undef __FUNCT__ 1512 #define __FUNCT__ "SNESDefaultRhsBC" 1513 /*@ 1514 SNESDefaultRhsBC - The default boundary condition function which does nothing. 1515 1516 Not collective 1517 1518 Input Parameters: 1519 . snes - The nonlinear solver context 1520 . rhs - The Rhs 1521 . ctx - The user-context 1522 1523 Level: beginner 1524 1525 .keywords: SNES, Rhs, boundary conditions 1526 .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate() 1527 @*/ 1528 PetscErrorCode SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx) 1529 { 1530 PetscFunctionBegin; 1531 PetscFunctionReturn(0); 1532 } 1533 1534 #undef __FUNCT__ 1535 #define __FUNCT__ "SNESSetSolutionBC" 1536 /*@C 1537 SNESSetSolutionBC - Sets the function which applies boundary conditions 1538 to the solution of each system. 1539 1540 Collective on SNES 1541 1542 Input Parameters: 1543 . snes - The nonlinear solver context 1544 . func - The function 1545 1546 Calling sequence of func: 1547 . func (SNES snes, Vec rsol, void *ctx); 1548 1549 . sol - The current solution vector 1550 . ctx - The user-context 1551 1552 Level: intermediate 1553 1554 .keywords: SNES, solution, boundary conditions 1555 .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate() 1556 @*/ 1557 PetscErrorCode SNESSetSolutionBC(SNES snes, PetscErrorCode (*func)(SNES, Vec, void *)) 1558 { 1559 PetscFunctionBegin; 1560 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1561 snes->applysolbc = func; 1562 PetscFunctionReturn(0); 1563 } 1564 1565 #undef __FUNCT__ 1566 #define __FUNCT__ "SNESDefaultSolutionBC" 1567 /*@ 1568 SNESDefaultSolutionBC - The default boundary condition function which does nothing. 1569 1570 Not collective 1571 1572 Input Parameters: 1573 . snes - The nonlinear solver context 1574 . sol - The solution 1575 . ctx - The user-context 1576 1577 Level: beginner 1578 1579 .keywords: SNES, solution, boundary conditions 1580 .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate() 1581 @*/ 1582 PetscErrorCode SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx) 1583 { 1584 PetscFunctionBegin; 1585 PetscFunctionReturn(0); 1586 } 1587 1588 #undef __FUNCT__ 1589 #define __FUNCT__ "SNESSetUpdate" 1590 /*@C 1591 SNESSetUpdate - Sets the general-purpose update function called 1592 at the beginning of every step of the iteration. 1593 1594 Collective on SNES 1595 1596 Input Parameters: 1597 . snes - The nonlinear solver context 1598 . func - The function 1599 1600 Calling sequence of func: 1601 . func (TS ts, int step); 1602 1603 . step - The current step of the iteration 1604 1605 Level: intermediate 1606 1607 .keywords: SNES, update 1608 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 1609 @*/ 1610 PetscErrorCode SNESSetUpdate(SNES snes, PetscErrorCode (*func)(SNES, PetscInt)) 1611 { 1612 PetscFunctionBegin; 1613 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1614 snes->update = func; 1615 PetscFunctionReturn(0); 1616 } 1617 1618 #undef __FUNCT__ 1619 #define __FUNCT__ "SNESDefaultUpdate" 1620 /*@ 1621 SNESDefaultUpdate - The default update function which does nothing. 1622 1623 Not collective 1624 1625 Input Parameters: 1626 . snes - The nonlinear solver context 1627 . step - The current step of the iteration 1628 1629 Level: intermediate 1630 1631 .keywords: SNES, update 1632 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 1633 @*/ 1634 PetscErrorCode SNESDefaultUpdate(SNES snes, PetscInt step) 1635 { 1636 PetscFunctionBegin; 1637 PetscFunctionReturn(0); 1638 } 1639 1640 #undef __FUNCT__ 1641 #define __FUNCT__ "SNESScaleStep_Private" 1642 /* 1643 SNESScaleStep_Private - Scales a step so that its length is less than the 1644 positive parameter delta. 1645 1646 Input Parameters: 1647 + snes - the SNES context 1648 . y - approximate solution of linear system 1649 . fnorm - 2-norm of current function 1650 - delta - trust region size 1651 1652 Output Parameters: 1653 + gpnorm - predicted function norm at the new point, assuming local 1654 linearization. The value is zero if the step lies within the trust 1655 region, and exceeds zero otherwise. 1656 - ynorm - 2-norm of the step 1657 1658 Note: 1659 For non-trust region methods such as SNESLS, the parameter delta 1660 is set to be the maximum allowable step size. 1661 1662 .keywords: SNES, nonlinear, scale, step 1663 */ 1664 PetscErrorCode SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 1665 { 1666 PetscReal nrm; 1667 PetscScalar cnorm; 1668 PetscErrorCode ierr; 1669 1670 PetscFunctionBegin; 1671 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1672 PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1673 PetscCheckSameComm(snes,1,y,2); 1674 1675 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1676 if (nrm > *delta) { 1677 nrm = *delta/nrm; 1678 *gpnorm = (1.0 - nrm)*(*fnorm); 1679 cnorm = nrm; 1680 ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 1681 *ynorm = *delta; 1682 } else { 1683 *gpnorm = 0.0; 1684 *ynorm = nrm; 1685 } 1686 PetscFunctionReturn(0); 1687 } 1688 1689 static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero", 1690 "not currently used", 1691 "line search failed", 1692 "reach maximum number of iterations", 1693 "function norm became NaN (not a number)", 1694 "not currently used", 1695 "number of function computations exceeded", 1696 "not currently used", 1697 "still iterating", 1698 "not currently used", 1699 "absolute size of function norm", 1700 "relative decrease in function norm", 1701 "step size is small", 1702 "not currently used", 1703 "not currently used", 1704 "small trust region"}; 1705 1706 #undef __FUNCT__ 1707 #define __FUNCT__ "SNESSolve" 1708 /*@ 1709 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 1710 SNESCreate() and optional routines of the form SNESSetXXX(). 1711 1712 Collective on SNES 1713 1714 Input Parameters: 1715 + snes - the SNES context 1716 - x - the solution vector 1717 1718 Notes: 1719 The user should initialize the vector,x, with the initial guess 1720 for the nonlinear solve prior to calling SNESSolve. In particular, 1721 to employ an initial guess of zero, the user should explicitly set 1722 this vector to zero by calling VecSet(). 1723 1724 Level: beginner 1725 1726 .keywords: SNES, nonlinear, solve 1727 1728 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs() 1729 @*/ 1730 PetscErrorCode SNESSolve(SNES snes,Vec x) 1731 { 1732 PetscErrorCode ierr; 1733 PetscTruth flg; 1734 1735 PetscFunctionBegin; 1736 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1737 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1738 PetscCheckSameComm(snes,1,x,2); 1739 if (!snes->solve) SETERRQ(PETSC_ERR_ORDER,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 1740 1741 if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1742 else {snes->vec_sol = snes->vec_sol_always = x;} 1743 if (snes->conv_hist_reset) snes->conv_hist_len = 0; 1744 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1745 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1746 ierr = (*(snes)->solve)(snes);CHKERRQ(ierr); 1747 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1748 ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 1749 if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1750 ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1751 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 1752 if (snes->printreason) { 1753 if (snes->reason > 0) { 1754 ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 1755 } else { 1756 ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 1757 } 1758 } 1759 1760 PetscFunctionReturn(0); 1761 } 1762 1763 /* --------- Internal routines for SNES Package --------- */ 1764 1765 #undef __FUNCT__ 1766 #define __FUNCT__ "SNESSetType" 1767 /*@C 1768 SNESSetType - Sets the method for the nonlinear solver. 1769 1770 Collective on SNES 1771 1772 Input Parameters: 1773 + snes - the SNES context 1774 - type - a known method 1775 1776 Options Database Key: 1777 . -snes_type <type> - Sets the method; use -help for a list 1778 of available methods (for instance, ls or tr) 1779 1780 Notes: 1781 See "petsc/include/petscsnes.h" for available methods (for instance) 1782 + SNESLS - Newton's method with line search 1783 (systems of nonlinear equations) 1784 . SNESTR - Newton's method with trust region 1785 (systems of nonlinear equations) 1786 1787 Normally, it is best to use the SNESSetFromOptions() command and then 1788 set the SNES solver type from the options database rather than by using 1789 this routine. Using the options database provides the user with 1790 maximum flexibility in evaluating the many nonlinear solvers. 1791 The SNESSetType() routine is provided for those situations where it 1792 is necessary to set the nonlinear solver independently of the command 1793 line or options database. This might be the case, for example, when 1794 the choice of solver changes during the execution of the program, 1795 and the user's application is taking responsibility for choosing the 1796 appropriate method. 1797 1798 Level: intermediate 1799 1800 .keywords: SNES, set, type 1801 1802 .seealso: SNESType, SNESCreate() 1803 1804 @*/ 1805 PetscErrorCode SNESSetType(SNES snes,const SNESType type) 1806 { 1807 PetscErrorCode ierr,(*r)(SNES); 1808 PetscTruth match; 1809 1810 PetscFunctionBegin; 1811 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1812 PetscValidCharPointer(type,2); 1813 1814 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 1815 if (match) PetscFunctionReturn(0); 1816 1817 if (snes->setupcalled) { 1818 ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 1819 snes->data = 0; 1820 } 1821 1822 /* Get the function pointers for the iterative method requested */ 1823 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1824 ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1825 if (!r) SETERRQ1(PETSC_ERR_ARG_UNKNOWN_TYPE,"Unable to find requested SNES type %s",type); 1826 if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 1827 snes->data = 0; 1828 ierr = (*r)(snes);CHKERRQ(ierr); 1829 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 1830 PetscFunctionReturn(0); 1831 } 1832 1833 1834 /* --------------------------------------------------------------------- */ 1835 #undef __FUNCT__ 1836 #define __FUNCT__ "SNESRegisterDestroy" 1837 /*@C 1838 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1839 registered by SNESRegisterDynamic(). 1840 1841 Not Collective 1842 1843 Level: advanced 1844 1845 .keywords: SNES, nonlinear, register, destroy 1846 1847 .seealso: SNESRegisterAll(), SNESRegisterAll() 1848 @*/ 1849 PetscErrorCode SNESRegisterDestroy(void) 1850 { 1851 PetscErrorCode ierr; 1852 1853 PetscFunctionBegin; 1854 if (SNESList) { 1855 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 1856 SNESList = 0; 1857 } 1858 SNESRegisterAllCalled = PETSC_FALSE; 1859 PetscFunctionReturn(0); 1860 } 1861 1862 #undef __FUNCT__ 1863 #define __FUNCT__ "SNESGetType" 1864 /*@C 1865 SNESGetType - Gets the SNES method type and name (as a string). 1866 1867 Not Collective 1868 1869 Input Parameter: 1870 . snes - nonlinear solver context 1871 1872 Output Parameter: 1873 . type - SNES method (a character string) 1874 1875 Level: intermediate 1876 1877 .keywords: SNES, nonlinear, get, type, name 1878 @*/ 1879 PetscErrorCode SNESGetType(SNES snes,SNESType *type) 1880 { 1881 PetscFunctionBegin; 1882 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1883 PetscValidPointer(type,2); 1884 *type = snes->type_name; 1885 PetscFunctionReturn(0); 1886 } 1887 1888 #undef __FUNCT__ 1889 #define __FUNCT__ "SNESGetSolution" 1890 /*@C 1891 SNESGetSolution - Returns the vector where the approximate solution is 1892 stored. 1893 1894 Not Collective, but Vec is parallel if SNES is parallel 1895 1896 Input Parameter: 1897 . snes - the SNES context 1898 1899 Output Parameter: 1900 . x - the solution 1901 1902 Level: advanced 1903 1904 .keywords: SNES, nonlinear, get, solution 1905 1906 .seealso: SNESGetFunction(), SNESGetSolutionUpdate() 1907 @*/ 1908 PetscErrorCode SNESGetSolution(SNES snes,Vec *x) 1909 { 1910 PetscFunctionBegin; 1911 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1912 PetscValidPointer(x,2); 1913 *x = snes->vec_sol_always; 1914 PetscFunctionReturn(0); 1915 } 1916 1917 #undef __FUNCT__ 1918 #define __FUNCT__ "SNESGetSolutionUpdate" 1919 /*@C 1920 SNESGetSolutionUpdate - Returns the vector where the solution update is 1921 stored. 1922 1923 Not Collective, but Vec is parallel if SNES is parallel 1924 1925 Input Parameter: 1926 . snes - the SNES context 1927 1928 Output Parameter: 1929 . x - the solution update 1930 1931 Level: advanced 1932 1933 .keywords: SNES, nonlinear, get, solution, update 1934 1935 .seealso: SNESGetSolution(), SNESGetFunction 1936 @*/ 1937 PetscErrorCode SNESGetSolutionUpdate(SNES snes,Vec *x) 1938 { 1939 PetscFunctionBegin; 1940 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1941 PetscValidPointer(x,2); 1942 *x = snes->vec_sol_update_always; 1943 PetscFunctionReturn(0); 1944 } 1945 1946 #undef __FUNCT__ 1947 #define __FUNCT__ "SNESGetFunction" 1948 /*@C 1949 SNESGetFunction - Returns the vector where the function is stored. 1950 1951 Not Collective, but Vec is parallel if SNES is parallel 1952 1953 Input Parameter: 1954 . snes - the SNES context 1955 1956 Output Parameter: 1957 + r - the function (or PETSC_NULL) 1958 . ctx - the function context (or PETSC_NULL) 1959 - func - the function (or PETSC_NULL) 1960 1961 Level: advanced 1962 1963 .keywords: SNES, nonlinear, get, function 1964 1965 .seealso: SNESSetFunction(), SNESGetSolution() 1966 @*/ 1967 PetscErrorCode SNESGetFunction(SNES snes,Vec *r,void **ctx,PetscErrorCode (**func)(SNES,Vec,Vec,void*)) 1968 { 1969 PetscFunctionBegin; 1970 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1971 if (r) *r = snes->vec_func_always; 1972 if (ctx) *ctx = snes->funP; 1973 if (func) *func = snes->computefunction; 1974 PetscFunctionReturn(0); 1975 } 1976 1977 #undef __FUNCT__ 1978 #define __FUNCT__ "SNESSetOptionsPrefix" 1979 /*@C 1980 SNESSetOptionsPrefix - Sets the prefix used for searching for all 1981 SNES options in the database. 1982 1983 Collective on SNES 1984 1985 Input Parameter: 1986 + snes - the SNES context 1987 - prefix - the prefix to prepend to all option names 1988 1989 Notes: 1990 A hyphen (-) must NOT be given at the beginning of the prefix name. 1991 The first character of all runtime options is AUTOMATICALLY the hyphen. 1992 1993 Level: advanced 1994 1995 .keywords: SNES, set, options, prefix, database 1996 1997 .seealso: SNESSetFromOptions() 1998 @*/ 1999 PetscErrorCode SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2000 { 2001 PetscErrorCode ierr; 2002 2003 PetscFunctionBegin; 2004 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2005 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2006 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2007 PetscFunctionReturn(0); 2008 } 2009 2010 #undef __FUNCT__ 2011 #define __FUNCT__ "SNESAppendOptionsPrefix" 2012 /*@C 2013 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2014 SNES options in the database. 2015 2016 Collective on SNES 2017 2018 Input Parameters: 2019 + snes - the SNES context 2020 - prefix - the prefix to prepend to all option names 2021 2022 Notes: 2023 A hyphen (-) must NOT be given at the beginning of the prefix name. 2024 The first character of all runtime options is AUTOMATICALLY the hyphen. 2025 2026 Level: advanced 2027 2028 .keywords: SNES, append, options, prefix, database 2029 2030 .seealso: SNESGetOptionsPrefix() 2031 @*/ 2032 PetscErrorCode SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2033 { 2034 PetscErrorCode ierr; 2035 2036 PetscFunctionBegin; 2037 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2038 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2039 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2040 PetscFunctionReturn(0); 2041 } 2042 2043 #undef __FUNCT__ 2044 #define __FUNCT__ "SNESGetOptionsPrefix" 2045 /*@C 2046 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2047 SNES options in the database. 2048 2049 Not Collective 2050 2051 Input Parameter: 2052 . snes - the SNES context 2053 2054 Output Parameter: 2055 . prefix - pointer to the prefix string used 2056 2057 Notes: On the fortran side, the user should pass in a string 'prifix' of 2058 sufficient length to hold the prefix. 2059 2060 Level: advanced 2061 2062 .keywords: SNES, get, options, prefix, database 2063 2064 .seealso: SNESAppendOptionsPrefix() 2065 @*/ 2066 PetscErrorCode SNESGetOptionsPrefix(SNES snes,char *prefix[]) 2067 { 2068 PetscErrorCode ierr; 2069 2070 PetscFunctionBegin; 2071 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2072 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2073 PetscFunctionReturn(0); 2074 } 2075 2076 2077 #undef __FUNCT__ 2078 #define __FUNCT__ "SNESRegister" 2079 /*@C 2080 SNESRegister - See SNESRegisterDynamic() 2081 2082 Level: advanced 2083 @*/ 2084 PetscErrorCode SNESRegister(const char sname[],const char path[],const char name[],PetscErrorCode (*function)(SNES)) 2085 { 2086 char fullname[PETSC_MAX_PATH_LEN]; 2087 PetscErrorCode ierr; 2088 2089 PetscFunctionBegin; 2090 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2091 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2092 PetscFunctionReturn(0); 2093 } 2094 2095 #undef __FUNCT__ 2096 #define __FUNCT__ "SNESTestLocalMin" 2097 PetscErrorCode SNESTestLocalMin(SNES snes) 2098 { 2099 PetscErrorCode ierr; 2100 PetscInt N,i,j; 2101 Vec u,uh,fh; 2102 PetscScalar value; 2103 PetscReal norm; 2104 2105 PetscFunctionBegin; 2106 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2107 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2108 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2109 2110 /* currently only works for sequential */ 2111 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2112 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2113 for (i=0; i<N; i++) { 2114 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2115 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %D\n",i);CHKERRQ(ierr); 2116 for (j=-10; j<11; j++) { 2117 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2118 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2119 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2120 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2121 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %D %18.16e\n",j,norm);CHKERRQ(ierr); 2122 value = -value; 2123 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2124 } 2125 } 2126 ierr = VecDestroy(uh);CHKERRQ(ierr); 2127 ierr = VecDestroy(fh);CHKERRQ(ierr); 2128 PetscFunctionReturn(0); 2129 } 2130