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