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