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