1 /*$Id: snes.c,v 1.235 2001/08/21 21:03:49 bsmith Exp $*/ 2 3 #include "src/snes/snesimpl.h" /*I "petscsnes.h" I*/ 4 5 PetscTruth SNESRegisterAllCalled = PETSC_FALSE; 6 PetscFList SNESList = PETSC_NULL; 7 8 /* Logging support */ 9 int SNES_COOKIE = 0; 10 int SNES_Solve = 0, SNES_LineSearch = 0, SNES_FunctionEval = 0, SNES_JacobianEval = 0; 11 12 #undef __FUNCT__ 13 #define __FUNCT__ "SNESView" 14 /*@C 15 SNESView - Prints the SNES data structure. 16 17 Collective on SNES 18 19 Input Parameters: 20 + SNES - the SNES context 21 - viewer - visualization context 22 23 Options Database Key: 24 . -snes_view - Calls SNESView() at end of SNESSolve() 25 26 Notes: 27 The available visualization contexts include 28 + PETSC_VIEWER_STDOUT_SELF - standard output (default) 29 - PETSC_VIEWER_STDOUT_WORLD - synchronized standard 30 output where only the first processor opens 31 the file. All other processors send their 32 data to the first processor to print. 33 34 The user can open an alternative visualization context with 35 PetscViewerASCIIOpen() - output to a specified file. 36 37 Level: beginner 38 39 .keywords: SNES, view 40 41 .seealso: PetscViewerASCIIOpen() 42 @*/ 43 int SNESView(SNES snes,PetscViewer viewer) 44 { 45 SNES_KSP_EW_ConvCtx *kctx; 46 int ierr; 47 KSP ksp; 48 char *type; 49 PetscTruth isascii,isstring; 50 51 PetscFunctionBegin; 52 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 53 if (!viewer) viewer = PETSC_VIEWER_STDOUT_(snes->comm); 54 PetscValidHeaderSpecific(viewer,PETSC_VIEWER_COOKIE,2); 55 PetscCheckSameComm(snes,1,viewer,2); 56 57 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);CHKERRQ(ierr); 58 ierr = PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_STRING,&isstring);CHKERRQ(ierr); 59 if (isascii) { 60 if (snes->prefix) { 61 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:(%s)\n",snes->prefix);CHKERRQ(ierr); 62 } else { 63 ierr = PetscViewerASCIIPrintf(viewer,"SNES Object:\n");CHKERRQ(ierr); 64 } 65 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 66 if (type) { 67 ierr = PetscViewerASCIIPrintf(viewer," type: %s\n",type);CHKERRQ(ierr); 68 } else { 69 ierr = PetscViewerASCIIPrintf(viewer," type: not set yet\n");CHKERRQ(ierr); 70 } 71 if (snes->view) { 72 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 73 ierr = (*snes->view)(snes,viewer);CHKERRQ(ierr); 74 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 75 } 76 ierr = PetscViewerASCIIPrintf(viewer," maximum iterations=%d, maximum function evaluations=%d\n",snes->max_its,snes->max_funcs);CHKERRQ(ierr); 77 ierr = PetscViewerASCIIPrintf(viewer," tolerances: relative=%g, absolute=%g, solution=%g\n", 78 snes->rtol,snes->atol,snes->xtol);CHKERRQ(ierr); 79 ierr = PetscViewerASCIIPrintf(viewer," total number of linear solver iterations=%d\n",snes->linear_its);CHKERRQ(ierr); 80 ierr = PetscViewerASCIIPrintf(viewer," total number of function evaluations=%d\n",snes->nfuncs);CHKERRQ(ierr); 81 if (snes->ksp_ewconv) { 82 kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 83 if (kctx) { 84 ierr = PetscViewerASCIIPrintf(viewer," Eisenstat-Walker computation of KSP relative tolerance (version %d)\n",kctx->version);CHKERRQ(ierr); 85 ierr = PetscViewerASCIIPrintf(viewer," rtol_0=%g, rtol_max=%g, threshold=%g\n",kctx->rtol_0,kctx->rtol_max,kctx->threshold);CHKERRQ(ierr); 86 ierr = PetscViewerASCIIPrintf(viewer," gamma=%g, alpha=%g, alpha2=%g\n",kctx->gamma,kctx->alpha,kctx->alpha2);CHKERRQ(ierr); 87 } 88 } 89 } else if (isstring) { 90 ierr = SNESGetType(snes,&type);CHKERRQ(ierr); 91 ierr = PetscViewerStringSPrintf(viewer," %-3.3s",type);CHKERRQ(ierr); 92 } 93 ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr); 94 ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 95 ierr = KSPView(ksp,viewer);CHKERRQ(ierr); 96 ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 97 PetscFunctionReturn(0); 98 } 99 100 /* 101 We retain a list of functions that also take SNES command 102 line options. These are called at the end SNESSetFromOptions() 103 */ 104 #define MAXSETFROMOPTIONS 5 105 static int numberofsetfromoptions; 106 static int (*othersetfromoptions[MAXSETFROMOPTIONS])(SNES); 107 108 #undef __FUNCT__ 109 #define __FUNCT__ "SNESAddOptionsChecker" 110 /*@C 111 SNESAddOptionsChecker - Adds an additional function to check for SNES options. 112 113 Not Collective 114 115 Input Parameter: 116 . snescheck - function that checks for options 117 118 Level: developer 119 120 .seealso: SNESSetFromOptions() 121 @*/ 122 int SNESAddOptionsChecker(int (*snescheck)(SNES)) 123 { 124 PetscFunctionBegin; 125 if (numberofsetfromoptions >= MAXSETFROMOPTIONS) { 126 SETERRQ1(PETSC_ERR_ARG_OUTOFRANGE, "Too many options checkers, only %d allowed", MAXSETFROMOPTIONS); 127 } 128 129 othersetfromoptions[numberofsetfromoptions++] = snescheck; 130 PetscFunctionReturn(0); 131 } 132 133 #undef __FUNCT__ 134 #define __FUNCT__ "SNESSetFromOptions" 135 /*@ 136 SNESSetFromOptions - Sets various SNES and KSP parameters from user options. 137 138 Collective on SNES 139 140 Input Parameter: 141 . snes - the SNES context 142 143 Options Database Keys: 144 + -snes_type <type> - ls, tr, umls, umtr, test 145 . -snes_stol - convergence tolerance in terms of the norm 146 of the change in the solution between steps 147 . -snes_atol <atol> - absolute tolerance of residual norm 148 . -snes_rtol <rtol> - relative decrease in tolerance norm from initial 149 . -snes_max_it <max_it> - maximum number of iterations 150 . -snes_max_funcs <max_funcs> - maximum number of function evaluations 151 . -snes_max_fail <max_fail> - maximum number of failures 152 . -snes_trtol <trtol> - trust region tolerance 153 . -snes_no_convergence_test - skip convergence test in nonlinear or minimization 154 solver; hence iterations will continue until max_it 155 or some other criterion is reached. Saves expense 156 of convergence test 157 . -snes_monitor - prints residual norm at each iteration 158 . -snes_vecmonitor - plots solution at each iteration 159 . -snes_vecmonitor_update - plots update to solution at each iteration 160 . -snes_xmonitor - plots residual norm at each iteration 161 . -snes_fd - use finite differences to compute Jacobian; very slow, only for testing 162 . -snes_mf_ksp_monitor - if using matrix-free multiply then print h at each KSP iteration 163 - -snes_print_converged_reason - print the reason for convergence/divergence after each solve 164 165 Options Database for Eisenstat-Walker method: 166 + -snes_ksp_ew_conv - use Eisenstat-Walker method for determining linear system convergence 167 . -snes_ksp_ew_version ver - version of Eisenstat-Walker method 168 . -snes_ksp_ew_rtol0 <rtol0> - Sets rtol0 169 . -snes_ksp_ew_rtolmax <rtolmax> - Sets rtolmax 170 . -snes_ksp_ew_gamma <gamma> - Sets gamma 171 . -snes_ksp_ew_alpha <alpha> - Sets alpha 172 . -snes_ksp_ew_alpha2 <alpha2> - Sets alpha2 173 - -snes_ksp_ew_threshold <threshold> - Sets threshold 174 175 Notes: 176 To see all options, run your program with the -help option or consult 177 the users manual. 178 179 Level: beginner 180 181 .keywords: SNES, nonlinear, set, options, database 182 183 .seealso: SNESSetOptionsPrefix() 184 @*/ 185 int SNESSetFromOptions(SNES snes) 186 { 187 KSP ksp; 188 SNES_KSP_EW_ConvCtx *kctx = (SNES_KSP_EW_ConvCtx *)snes->kspconvctx; 189 PetscTruth flg; 190 int ierr, 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 int 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 int 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 int 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 int 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 int 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 int 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 int 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 int 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 int 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 int SNESPublish_Petsc(PetscObject obj) 548 { 549 #if defined(PETSC_HAVE_AMS) 550 SNES v = (SNES) obj; 551 int 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 int SNESCreate(MPI_Comm comm,SNES *outsnes) 599 { 600 int 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 int SNESSetFunction(SNES snes,Vec r,int (*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 int SNESSetRhs(SNES snes,Vec rhs) 730 { 731 int 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 int SNESComputeFunction(SNES snes,Vec x,Vec y) 774 { 775 int 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 int SNESComputeJacobian(SNES snes,Vec X,Mat *A,Mat *B,MatStructure *flg) 828 { 829 int 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 int SNESSetJacobian(SNES snes,Mat A,Mat B,int (*func)(SNES,Vec,Mat*,Mat*,MatStructure*,void*),void *ctx) 892 { 893 int 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 int SNESGetJacobian(SNES snes,Mat *A,Mat *B,void **ctx,int (**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 int 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 int SNESSetUp(SNES snes,Vec x) 977 { 978 int 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 int SNESDestroy(SNES snes) 1074 { 1075 int i,ierr; 1076 1077 PetscFunctionBegin; 1078 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1079 if (--snes->refct > 0) PetscFunctionReturn(0); 1080 1081 /* if memory was published with AMS then destroy it */ 1082 ierr = PetscObjectDepublish(snes);CHKERRQ(ierr); 1083 1084 if (snes->destroy) {ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr);} 1085 if (snes->kspconvctx) {ierr = PetscFree(snes->kspconvctx);CHKERRQ(ierr);} 1086 if (snes->jacobian) {ierr = MatDestroy(snes->jacobian);CHKERRQ(ierr);} 1087 if (snes->jacobian_pre) {ierr = MatDestroy(snes->jacobian_pre);CHKERRQ(ierr);} 1088 if (snes->afine) {ierr = VecDestroy(snes->afine);CHKERRQ(ierr);} 1089 ierr = KSPDestroy(snes->ksp);CHKERRQ(ierr); 1090 if (snes->vwork) {ierr = VecDestroyVecs(snes->vwork,snes->nvwork);CHKERRQ(ierr);} 1091 for (i=0; i<snes->numbermonitors; i++) { 1092 if (snes->monitordestroy[i]) { 1093 ierr = (*snes->monitordestroy[i])(snes->monitorcontext[i]);CHKERRQ(ierr); 1094 } 1095 } 1096 PetscLogObjectDestroy((PetscObject)snes); 1097 PetscHeaderDestroy((PetscObject)snes); 1098 PetscFunctionReturn(0); 1099 } 1100 1101 /* ----------- Routines to set solver parameters ---------- */ 1102 1103 #undef __FUNCT__ 1104 #define __FUNCT__ "SNESSetTolerances" 1105 /*@ 1106 SNESSetTolerances - Sets various parameters used in convergence tests. 1107 1108 Collective on SNES 1109 1110 Input Parameters: 1111 + snes - the SNES context 1112 . atol - absolute convergence tolerance 1113 . rtol - relative convergence tolerance 1114 . stol - convergence tolerance in terms of the norm 1115 of the change in the solution between steps 1116 . maxit - maximum number of iterations 1117 - maxf - maximum number of function evaluations 1118 1119 Options Database Keys: 1120 + -snes_atol <atol> - Sets atol 1121 . -snes_rtol <rtol> - Sets rtol 1122 . -snes_stol <stol> - Sets stol 1123 . -snes_max_it <maxit> - Sets maxit 1124 - -snes_max_funcs <maxf> - Sets maxf 1125 1126 Notes: 1127 The default maximum number of iterations is 50. 1128 The default maximum number of function evaluations is 1000. 1129 1130 Level: intermediate 1131 1132 .keywords: SNES, nonlinear, set, convergence, tolerances 1133 1134 .seealso: SNESSetTrustRegionTolerance(), SNESSetMinimizationFunctionTolerance() 1135 @*/ 1136 int SNESSetTolerances(SNES snes,PetscReal atol,PetscReal rtol,PetscReal stol,int maxit,int maxf) 1137 { 1138 PetscFunctionBegin; 1139 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1140 if (atol != PETSC_DEFAULT) snes->atol = atol; 1141 if (rtol != PETSC_DEFAULT) snes->rtol = rtol; 1142 if (stol != PETSC_DEFAULT) snes->xtol = stol; 1143 if (maxit != PETSC_DEFAULT) snes->max_its = maxit; 1144 if (maxf != PETSC_DEFAULT) snes->max_funcs = maxf; 1145 PetscFunctionReturn(0); 1146 } 1147 1148 #undef __FUNCT__ 1149 #define __FUNCT__ "SNESGetTolerances" 1150 /*@ 1151 SNESGetTolerances - Gets various parameters used in convergence tests. 1152 1153 Not Collective 1154 1155 Input Parameters: 1156 + snes - the SNES context 1157 . atol - absolute convergence tolerance 1158 . rtol - relative convergence tolerance 1159 . stol - convergence tolerance in terms of the norm 1160 of the change in the solution between steps 1161 . maxit - maximum number of iterations 1162 - maxf - maximum number of function evaluations 1163 1164 Notes: 1165 The user can specify PETSC_NULL for any parameter that is not needed. 1166 1167 Level: intermediate 1168 1169 .keywords: SNES, nonlinear, get, convergence, tolerances 1170 1171 .seealso: SNESSetTolerances() 1172 @*/ 1173 int SNESGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,int *maxit,int *maxf) 1174 { 1175 PetscFunctionBegin; 1176 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1177 if (atol) *atol = snes->atol; 1178 if (rtol) *rtol = snes->rtol; 1179 if (stol) *stol = snes->xtol; 1180 if (maxit) *maxit = snes->max_its; 1181 if (maxf) *maxf = snes->max_funcs; 1182 PetscFunctionReturn(0); 1183 } 1184 1185 #undef __FUNCT__ 1186 #define __FUNCT__ "SNESSetTrustRegionTolerance" 1187 /*@ 1188 SNESSetTrustRegionTolerance - Sets the trust region parameter tolerance. 1189 1190 Collective on SNES 1191 1192 Input Parameters: 1193 + snes - the SNES context 1194 - tol - tolerance 1195 1196 Options Database Key: 1197 . -snes_trtol <tol> - Sets tol 1198 1199 Level: intermediate 1200 1201 .keywords: SNES, nonlinear, set, trust region, tolerance 1202 1203 .seealso: SNESSetTolerances(), SNESSetMinimizationFunctionTolerance() 1204 @*/ 1205 int SNESSetTrustRegionTolerance(SNES snes,PetscReal tol) 1206 { 1207 PetscFunctionBegin; 1208 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1209 snes->deltatol = tol; 1210 PetscFunctionReturn(0); 1211 } 1212 1213 /* 1214 Duplicate the lg monitors for SNES from KSP; for some reason with 1215 dynamic libraries things don't work under Sun4 if we just use 1216 macros instead of functions 1217 */ 1218 #undef __FUNCT__ 1219 #define __FUNCT__ "SNESLGMonitor" 1220 int SNESLGMonitor(SNES snes,int it,PetscReal norm,void *ctx) 1221 { 1222 int ierr; 1223 1224 PetscFunctionBegin; 1225 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1226 ierr = KSPLGMonitor((KSP)snes,it,norm,ctx);CHKERRQ(ierr); 1227 PetscFunctionReturn(0); 1228 } 1229 1230 #undef __FUNCT__ 1231 #define __FUNCT__ "SNESLGMonitorCreate" 1232 int SNESLGMonitorCreate(const char host[],const char label[],int x,int y,int m,int n,PetscDrawLG *draw) 1233 { 1234 int ierr; 1235 1236 PetscFunctionBegin; 1237 ierr = KSPLGMonitorCreate(host,label,x,y,m,n,draw);CHKERRQ(ierr); 1238 PetscFunctionReturn(0); 1239 } 1240 1241 #undef __FUNCT__ 1242 #define __FUNCT__ "SNESLGMonitorDestroy" 1243 int SNESLGMonitorDestroy(PetscDrawLG draw) 1244 { 1245 int ierr; 1246 1247 PetscFunctionBegin; 1248 ierr = KSPLGMonitorDestroy(draw);CHKERRQ(ierr); 1249 PetscFunctionReturn(0); 1250 } 1251 1252 /* ------------ Routines to set performance monitoring options ----------- */ 1253 1254 #undef __FUNCT__ 1255 #define __FUNCT__ "SNESSetMonitor" 1256 /*@C 1257 SNESSetMonitor - Sets an ADDITIONAL function that is to be used at every 1258 iteration of the nonlinear solver to display the iteration's 1259 progress. 1260 1261 Collective on SNES 1262 1263 Input Parameters: 1264 + snes - the SNES context 1265 . func - monitoring routine 1266 . mctx - [optional] user-defined context for private data for the 1267 monitor routine (use PETSC_NULL if no context is desitre) 1268 - monitordestroy - [optional] routine that frees monitor context 1269 (may be PETSC_NULL) 1270 1271 Calling sequence of func: 1272 $ int func(SNES snes,int its, PetscReal norm,void *mctx) 1273 1274 + snes - the SNES context 1275 . its - iteration number 1276 . norm - 2-norm function value (may be estimated) 1277 - mctx - [optional] monitoring context 1278 1279 Options Database Keys: 1280 + -snes_monitor - sets SNESDefaultMonitor() 1281 . -snes_xmonitor - sets line graph monitor, 1282 uses SNESLGMonitorCreate() 1283 _ -snes_cancelmonitors - cancels all monitors that have 1284 been hardwired into a code by 1285 calls to SNESSetMonitor(), but 1286 does not cancel those set via 1287 the options database. 1288 1289 Notes: 1290 Several different monitoring routines may be set by calling 1291 SNESSetMonitor() multiple times; all will be called in the 1292 order in which they were set. 1293 1294 Level: intermediate 1295 1296 .keywords: SNES, nonlinear, set, monitor 1297 1298 .seealso: SNESDefaultMonitor(), SNESClearMonitor() 1299 @*/ 1300 int SNESSetMonitor(SNES snes,int (*func)(SNES,int,PetscReal,void*),void *mctx,int (*monitordestroy)(void *)) 1301 { 1302 PetscFunctionBegin; 1303 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1304 if (snes->numbermonitors >= MAXSNESMONITORS) { 1305 SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many monitors set"); 1306 } 1307 1308 snes->monitor[snes->numbermonitors] = func; 1309 snes->monitordestroy[snes->numbermonitors] = monitordestroy; 1310 snes->monitorcontext[snes->numbermonitors++] = (void*)mctx; 1311 PetscFunctionReturn(0); 1312 } 1313 1314 #undef __FUNCT__ 1315 #define __FUNCT__ "SNESClearMonitor" 1316 /*@C 1317 SNESClearMonitor - Clears all the monitor functions for a SNES object. 1318 1319 Collective on SNES 1320 1321 Input Parameters: 1322 . snes - the SNES context 1323 1324 Options Database Key: 1325 . -snes_cancelmonitors - cancels all monitors that have been hardwired 1326 into a code by calls to SNESSetMonitor(), but does not cancel those 1327 set via the options database 1328 1329 Notes: 1330 There is no way to clear one specific monitor from a SNES object. 1331 1332 Level: intermediate 1333 1334 .keywords: SNES, nonlinear, set, monitor 1335 1336 .seealso: SNESDefaultMonitor(), SNESSetMonitor() 1337 @*/ 1338 int SNESClearMonitor(SNES snes) 1339 { 1340 PetscFunctionBegin; 1341 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1342 snes->numbermonitors = 0; 1343 PetscFunctionReturn(0); 1344 } 1345 1346 #undef __FUNCT__ 1347 #define __FUNCT__ "SNESSetConvergenceTest" 1348 /*@C 1349 SNESSetConvergenceTest - Sets the function that is to be used 1350 to test for convergence of the nonlinear iterative solution. 1351 1352 Collective on SNES 1353 1354 Input Parameters: 1355 + snes - the SNES context 1356 . func - routine to test for convergence 1357 - cctx - [optional] context for private data for the convergence routine 1358 (may be PETSC_NULL) 1359 1360 Calling sequence of func: 1361 $ int func (SNES snes,PetscReal xnorm,PetscReal gnorm,PetscReal f,SNESConvergedReason *reason,void *cctx) 1362 1363 + snes - the SNES context 1364 . cctx - [optional] convergence context 1365 . reason - reason for convergence/divergence 1366 . xnorm - 2-norm of current iterate 1367 . gnorm - 2-norm of current step 1368 - f - 2-norm of function 1369 1370 Level: advanced 1371 1372 .keywords: SNES, nonlinear, set, convergence, test 1373 1374 .seealso: SNESConverged_LS(), SNESConverged_TR() 1375 @*/ 1376 int SNESSetConvergenceTest(SNES snes,int (*func)(SNES,PetscReal,PetscReal,PetscReal,SNESConvergedReason*,void*),void *cctx) 1377 { 1378 PetscFunctionBegin; 1379 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1380 (snes)->converged = func; 1381 (snes)->cnvP = cctx; 1382 PetscFunctionReturn(0); 1383 } 1384 1385 #undef __FUNCT__ 1386 #define __FUNCT__ "SNESGetConvergedReason" 1387 /*@C 1388 SNESGetConvergedReason - Gets the reason the SNES iteration was stopped. 1389 1390 Not Collective 1391 1392 Input Parameter: 1393 . snes - the SNES context 1394 1395 Output Parameter: 1396 . reason - negative value indicates diverged, positive value converged, see petscsnes.h or the 1397 manual pages for the individual convergence tests for complete lists 1398 1399 Level: intermediate 1400 1401 Notes: Can only be called after the call the SNESSolve() is complete. 1402 1403 .keywords: SNES, nonlinear, set, convergence, test 1404 1405 .seealso: SNESSetConvergenceTest(), SNESConverged_LS(), SNESConverged_TR(), SNESConvergedReason 1406 @*/ 1407 int SNESGetConvergedReason(SNES snes,SNESConvergedReason *reason) 1408 { 1409 PetscFunctionBegin; 1410 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1411 PetscValidPointer(reason,2); 1412 *reason = snes->reason; 1413 PetscFunctionReturn(0); 1414 } 1415 1416 #undef __FUNCT__ 1417 #define __FUNCT__ "SNESSetConvergenceHistory" 1418 /*@ 1419 SNESSetConvergenceHistory - Sets the array used to hold the convergence history. 1420 1421 Collective on SNES 1422 1423 Input Parameters: 1424 + snes - iterative context obtained from SNESCreate() 1425 . a - array to hold history 1426 . its - integer array holds the number of linear iterations for each solve. 1427 . na - size of a and its 1428 - reset - PETSC_TRUE indicates each new nonlinear solve resets the history counter to zero, 1429 else it continues storing new values for new nonlinear solves after the old ones 1430 1431 Notes: 1432 If set, this array will contain the function norms computed 1433 at each step. 1434 1435 This routine is useful, e.g., when running a code for purposes 1436 of accurate performance monitoring, when no I/O should be done 1437 during the section of code that is being timed. 1438 1439 Level: intermediate 1440 1441 .keywords: SNES, set, convergence, history 1442 1443 .seealso: SNESGetConvergenceHistory() 1444 1445 @*/ 1446 int SNESSetConvergenceHistory(SNES snes,PetscReal a[],int *its,int na,PetscTruth reset) 1447 { 1448 PetscFunctionBegin; 1449 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1450 if (na) PetscValidScalarPointer(a,2); 1451 snes->conv_hist = a; 1452 snes->conv_hist_its = its; 1453 snes->conv_hist_max = na; 1454 snes->conv_hist_reset = reset; 1455 PetscFunctionReturn(0); 1456 } 1457 1458 #undef __FUNCT__ 1459 #define __FUNCT__ "SNESGetConvergenceHistory" 1460 /*@C 1461 SNESGetConvergenceHistory - Gets the array used to hold the convergence history. 1462 1463 Collective on SNES 1464 1465 Input Parameter: 1466 . snes - iterative context obtained from SNESCreate() 1467 1468 Output Parameters: 1469 . a - array to hold history 1470 . its - integer array holds the number of linear iterations (or 1471 negative if not converged) for each solve. 1472 - na - size of a and its 1473 1474 Notes: 1475 The calling sequence for this routine in Fortran is 1476 $ call SNESGetConvergenceHistory(SNES snes, integer na, integer ierr) 1477 1478 This routine is useful, e.g., when running a code for purposes 1479 of accurate performance monitoring, when no I/O should be done 1480 during the section of code that is being timed. 1481 1482 Level: intermediate 1483 1484 .keywords: SNES, get, convergence, history 1485 1486 .seealso: SNESSetConvergencHistory() 1487 1488 @*/ 1489 int SNESGetConvergenceHistory(SNES snes,PetscReal *a[],int *its[],int *na) 1490 { 1491 PetscFunctionBegin; 1492 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1493 if (a) *a = snes->conv_hist; 1494 if (its) *its = snes->conv_hist_its; 1495 if (na) *na = snes->conv_hist_len; 1496 PetscFunctionReturn(0); 1497 } 1498 1499 #undef __FUNCT__ 1500 #define __FUNCT__ "SNESSetRhsBC" 1501 /*@C 1502 SNESSetRhsBC - Sets the function which applies boundary conditions 1503 to the Rhs of each system. 1504 1505 Collective on SNES 1506 1507 Input Parameters: 1508 . snes - The nonlinear solver context 1509 . func - The function 1510 1511 Calling sequence of func: 1512 . func (SNES snes, Vec rhs, void *ctx); 1513 1514 . rhs - The current rhs vector 1515 . ctx - The user-context 1516 1517 Level: intermediate 1518 1519 .keywords: SNES, Rhs, boundary conditions 1520 .seealso SNESDefaultRhsBC(), SNESSetSolutionBC(), SNESSetUpdate() 1521 @*/ 1522 int SNESSetRhsBC(SNES snes, int (*func)(SNES, Vec, void *)) 1523 { 1524 PetscFunctionBegin; 1525 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1526 snes->applyrhsbc = func; 1527 PetscFunctionReturn(0); 1528 } 1529 1530 #undef __FUNCT__ 1531 #define __FUNCT__ "SNESDefaultRhsBC" 1532 /*@ 1533 SNESDefaultRhsBC - The default boundary condition function which does nothing. 1534 1535 Not collective 1536 1537 Input Parameters: 1538 . snes - The nonlinear solver context 1539 . rhs - The Rhs 1540 . ctx - The user-context 1541 1542 Level: beginner 1543 1544 .keywords: SNES, Rhs, boundary conditions 1545 .seealso SNESSetRhsBC(), SNESDefaultSolutionBC(), SNESDefaultUpdate() 1546 @*/ 1547 int SNESDefaultRhsBC(SNES snes, Vec rhs, void *ctx) 1548 { 1549 PetscFunctionBegin; 1550 PetscFunctionReturn(0); 1551 } 1552 1553 #undef __FUNCT__ 1554 #define __FUNCT__ "SNESSetSolutionBC" 1555 /*@C 1556 SNESSetSolutionBC - Sets the function which applies boundary conditions 1557 to the solution of each system. 1558 1559 Collective on SNES 1560 1561 Input Parameters: 1562 . snes - The nonlinear solver context 1563 . func - The function 1564 1565 Calling sequence of func: 1566 . func (SNES snes, Vec rsol, void *ctx); 1567 1568 . sol - The current solution vector 1569 . ctx - The user-context 1570 1571 Level: intermediate 1572 1573 .keywords: SNES, solution, boundary conditions 1574 .seealso SNESDefaultSolutionBC(), SNESSetRhsBC(), SNESSetUpdate() 1575 @*/ 1576 int SNESSetSolutionBC(SNES snes, int (*func)(SNES, Vec, void *)) 1577 { 1578 PetscFunctionBegin; 1579 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1580 snes->applysolbc = func; 1581 PetscFunctionReturn(0); 1582 } 1583 1584 #undef __FUNCT__ 1585 #define __FUNCT__ "SNESDefaultSolutionBC" 1586 /*@ 1587 SNESDefaultSolutionBC - The default boundary condition function which does nothing. 1588 1589 Not collective 1590 1591 Input Parameters: 1592 . snes - The nonlinear solver context 1593 . sol - The solution 1594 . ctx - The user-context 1595 1596 Level: beginner 1597 1598 .keywords: SNES, solution, boundary conditions 1599 .seealso SNESSetSolutionBC(), SNESDefaultRhsBC(), SNESDefaultUpdate() 1600 @*/ 1601 int SNESDefaultSolutionBC(SNES snes, Vec sol, void *ctx) 1602 { 1603 PetscFunctionBegin; 1604 PetscFunctionReturn(0); 1605 } 1606 1607 #undef __FUNCT__ 1608 #define __FUNCT__ "SNESSetUpdate" 1609 /*@C 1610 SNESSetUpdate - Sets the general-purpose update function called 1611 at the beginning of every step of the iteration. 1612 1613 Collective on SNES 1614 1615 Input Parameters: 1616 . snes - The nonlinear solver context 1617 . func - The function 1618 1619 Calling sequence of func: 1620 . func (TS ts, int step); 1621 1622 . step - The current step of the iteration 1623 1624 Level: intermediate 1625 1626 .keywords: SNES, update 1627 .seealso SNESDefaultUpdate(), SNESSetRhsBC(), SNESSetSolutionBC() 1628 @*/ 1629 int SNESSetUpdate(SNES snes, int (*func)(SNES, int)) 1630 { 1631 PetscFunctionBegin; 1632 PetscValidHeaderSpecific(snes, SNES_COOKIE,1); 1633 snes->update = func; 1634 PetscFunctionReturn(0); 1635 } 1636 1637 #undef __FUNCT__ 1638 #define __FUNCT__ "SNESDefaultUpdate" 1639 /*@ 1640 SNESDefaultUpdate - The default update function which does nothing. 1641 1642 Not collective 1643 1644 Input Parameters: 1645 . snes - The nonlinear solver context 1646 . step - The current step of the iteration 1647 1648 Level: intermediate 1649 1650 .keywords: SNES, update 1651 .seealso SNESSetUpdate(), SNESDefaultRhsBC(), SNESDefaultSolutionBC() 1652 @*/ 1653 int SNESDefaultUpdate(SNES snes, int step) 1654 { 1655 PetscFunctionBegin; 1656 PetscFunctionReturn(0); 1657 } 1658 1659 #undef __FUNCT__ 1660 #define __FUNCT__ "SNESScaleStep_Private" 1661 /* 1662 SNESScaleStep_Private - Scales a step so that its length is less than the 1663 positive parameter delta. 1664 1665 Input Parameters: 1666 + snes - the SNES context 1667 . y - approximate solution of linear system 1668 . fnorm - 2-norm of current function 1669 - delta - trust region size 1670 1671 Output Parameters: 1672 + gpnorm - predicted function norm at the new point, assuming local 1673 linearization. The value is zero if the step lies within the trust 1674 region, and exceeds zero otherwise. 1675 - ynorm - 2-norm of the step 1676 1677 Note: 1678 For non-trust region methods such as SNESLS, the parameter delta 1679 is set to be the maximum allowable step size. 1680 1681 .keywords: SNES, nonlinear, scale, step 1682 */ 1683 int SNESScaleStep_Private(SNES snes,Vec y,PetscReal *fnorm,PetscReal *delta,PetscReal *gpnorm,PetscReal *ynorm) 1684 { 1685 PetscReal nrm; 1686 PetscScalar cnorm; 1687 int ierr; 1688 1689 PetscFunctionBegin; 1690 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1691 PetscValidHeaderSpecific(y,VEC_COOKIE,2); 1692 PetscCheckSameComm(snes,1,y,2); 1693 1694 ierr = VecNorm(y,NORM_2,&nrm);CHKERRQ(ierr); 1695 if (nrm > *delta) { 1696 nrm = *delta/nrm; 1697 *gpnorm = (1.0 - nrm)*(*fnorm); 1698 cnorm = nrm; 1699 ierr = VecScale(&cnorm,y);CHKERRQ(ierr); 1700 *ynorm = *delta; 1701 } else { 1702 *gpnorm = 0.0; 1703 *ynorm = nrm; 1704 } 1705 PetscFunctionReturn(0); 1706 } 1707 1708 static const char *convergedreasons[] = {"appears to located a local minimum instead of a zero", 1709 "not currently used", 1710 "line search failed", 1711 "reach maximum number of iterations", 1712 "function norm became NaN (not a number)", 1713 "not currently used", 1714 "number of function computations exceeded", 1715 "not currently used", 1716 "still iterating", 1717 "not currently used", 1718 "absolute size of function norm", 1719 "relative decrease in function norm", 1720 "step size is small", 1721 "not currently used", 1722 "not currently used", 1723 "small trust region"}; 1724 1725 #undef __FUNCT__ 1726 #define __FUNCT__ "SNESSolve" 1727 /*@ 1728 SNESSolve - Solves a nonlinear system. Call SNESSolve after calling 1729 SNESCreate() and optional routines of the form SNESSetXXX(). 1730 1731 Collective on SNES 1732 1733 Input Parameters: 1734 + snes - the SNES context 1735 - x - the solution vector 1736 1737 Notes: 1738 The user should initialize the vector,x, with the initial guess 1739 for the nonlinear solve prior to calling SNESSolve. In particular, 1740 to employ an initial guess of zero, the user should explicitly set 1741 this vector to zero by calling VecSet(). 1742 1743 Level: beginner 1744 1745 .keywords: SNES, nonlinear, solve 1746 1747 .seealso: SNESCreate(), SNESDestroy(), SNESSetFunction(), SNESSetJacobian(), SNESSetRhs() 1748 @*/ 1749 int SNESSolve(SNES snes,Vec x) 1750 { 1751 int ierr; 1752 PetscTruth flg; 1753 1754 PetscFunctionBegin; 1755 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1756 PetscValidHeaderSpecific(x,VEC_COOKIE,2); 1757 PetscCheckSameComm(snes,1,x,2); 1758 if (!snes->solve) SETERRQ(1,"SNESSetType() or SNESSetFromOptions() must be called before SNESSolve()"); 1759 1760 if (!snes->setupcalled) {ierr = SNESSetUp(snes,x);CHKERRQ(ierr);} 1761 else {snes->vec_sol = snes->vec_sol_always = x;} 1762 if (snes->conv_hist_reset == PETSC_TRUE) snes->conv_hist_len = 0; 1763 ierr = PetscLogEventBegin(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1764 snes->nfuncs = 0; snes->linear_its = 0; snes->numFailures = 0; 1765 ierr = (*(snes)->solve)(snes);CHKERRQ(ierr); 1766 ierr = PetscLogEventEnd(SNES_Solve,snes,0,0,0);CHKERRQ(ierr); 1767 ierr = PetscOptionsHasName(snes->prefix,"-snes_view",&flg);CHKERRQ(ierr); 1768 if (flg && !PetscPreLoadingOn) { ierr = SNESView(snes,PETSC_VIEWER_STDOUT_(snes->comm));CHKERRQ(ierr); } 1769 ierr = PetscOptionsHasName(snes->prefix,"-snes_test_local_min",&flg);CHKERRQ(ierr); 1770 if (flg && !PetscPreLoadingOn) { ierr = SNESTestLocalMin(snes);CHKERRQ(ierr); } 1771 if (snes->printreason) { 1772 if (snes->reason > 0) { 1773 ierr = PetscPrintf(snes->comm,"Nonlinear solve converged due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 1774 } else { 1775 ierr = PetscPrintf(snes->comm,"Nonlinear solve did not converge due to %s\n",convergedreasons[snes->reason+8]);CHKERRQ(ierr); 1776 } 1777 } 1778 1779 PetscFunctionReturn(0); 1780 } 1781 1782 /* --------- Internal routines for SNES Package --------- */ 1783 1784 #undef __FUNCT__ 1785 #define __FUNCT__ "SNESSetType" 1786 /*@C 1787 SNESSetType - Sets the method for the nonlinear solver. 1788 1789 Collective on SNES 1790 1791 Input Parameters: 1792 + snes - the SNES context 1793 - type - a known method 1794 1795 Options Database Key: 1796 . -snes_type <type> - Sets the method; use -help for a list 1797 of available methods (for instance, ls or tr) 1798 1799 Notes: 1800 See "petsc/include/petscsnes.h" for available methods (for instance) 1801 + SNESLS - Newton's method with line search 1802 (systems of nonlinear equations) 1803 . SNESTR - Newton's method with trust region 1804 (systems of nonlinear equations) 1805 1806 Normally, it is best to use the SNESSetFromOptions() command and then 1807 set the SNES solver type from the options database rather than by using 1808 this routine. Using the options database provides the user with 1809 maximum flexibility in evaluating the many nonlinear solvers. 1810 The SNESSetType() routine is provided for those situations where it 1811 is necessary to set the nonlinear solver independently of the command 1812 line or options database. This might be the case, for example, when 1813 the choice of solver changes during the execution of the program, 1814 and the user's application is taking responsibility for choosing the 1815 appropriate method. 1816 1817 Level: intermediate 1818 1819 .keywords: SNES, set, type 1820 1821 .seealso: SNESType, SNESCreate() 1822 1823 @*/ 1824 int SNESSetType(SNES snes,const SNESType type) 1825 { 1826 int ierr,(*r)(SNES); 1827 PetscTruth match; 1828 1829 PetscFunctionBegin; 1830 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1831 PetscValidCharPointer(type,2); 1832 1833 ierr = PetscTypeCompare((PetscObject)snes,type,&match);CHKERRQ(ierr); 1834 if (match) PetscFunctionReturn(0); 1835 1836 if (snes->setupcalled) { 1837 ierr = (*(snes)->destroy)(snes);CHKERRQ(ierr); 1838 snes->data = 0; 1839 } 1840 1841 /* Get the function pointers for the iterative method requested */ 1842 if (!SNESRegisterAllCalled) {ierr = SNESRegisterAll(PETSC_NULL);CHKERRQ(ierr);} 1843 1844 ierr = PetscFListFind(snes->comm,SNESList,type,(void (**)(void)) &r);CHKERRQ(ierr); 1845 1846 if (!r) SETERRQ1(1,"Unable to find requested SNES type %s",type); 1847 1848 if (snes->data) {ierr = PetscFree(snes->data);CHKERRQ(ierr);} 1849 snes->data = 0; 1850 ierr = (*r)(snes);CHKERRQ(ierr); 1851 ierr = PetscObjectChangeTypeName((PetscObject)snes,type);CHKERRQ(ierr); 1852 1853 PetscFunctionReturn(0); 1854 } 1855 1856 1857 /* --------------------------------------------------------------------- */ 1858 #undef __FUNCT__ 1859 #define __FUNCT__ "SNESRegisterDestroy" 1860 /*@C 1861 SNESRegisterDestroy - Frees the list of nonlinear solvers that were 1862 registered by SNESRegisterDynamic(). 1863 1864 Not Collective 1865 1866 Level: advanced 1867 1868 .keywords: SNES, nonlinear, register, destroy 1869 1870 .seealso: SNESRegisterAll(), SNESRegisterAll() 1871 @*/ 1872 int SNESRegisterDestroy(void) 1873 { 1874 int ierr; 1875 1876 PetscFunctionBegin; 1877 if (SNESList) { 1878 ierr = PetscFListDestroy(&SNESList);CHKERRQ(ierr); 1879 SNESList = 0; 1880 } 1881 SNESRegisterAllCalled = PETSC_FALSE; 1882 PetscFunctionReturn(0); 1883 } 1884 1885 #undef __FUNCT__ 1886 #define __FUNCT__ "SNESGetType" 1887 /*@C 1888 SNESGetType - Gets the SNES method type and name (as a string). 1889 1890 Not Collective 1891 1892 Input Parameter: 1893 . snes - nonlinear solver context 1894 1895 Output Parameter: 1896 . type - SNES method (a character string) 1897 1898 Level: intermediate 1899 1900 .keywords: SNES, nonlinear, get, type, name 1901 @*/ 1902 int SNESGetType(SNES snes,SNESType *type) 1903 { 1904 PetscFunctionBegin; 1905 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1906 PetscValidPointer(type,2); 1907 *type = snes->type_name; 1908 PetscFunctionReturn(0); 1909 } 1910 1911 #undef __FUNCT__ 1912 #define __FUNCT__ "SNESGetSolution" 1913 /*@C 1914 SNESGetSolution - Returns the vector where the approximate solution is 1915 stored. 1916 1917 Not Collective, but Vec is parallel if SNES is parallel 1918 1919 Input Parameter: 1920 . snes - the SNES context 1921 1922 Output Parameter: 1923 . x - the solution 1924 1925 Level: advanced 1926 1927 .keywords: SNES, nonlinear, get, solution 1928 1929 .seealso: SNESGetFunction(), SNESGetSolutionUpdate() 1930 @*/ 1931 int SNESGetSolution(SNES snes,Vec *x) 1932 { 1933 PetscFunctionBegin; 1934 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1935 PetscValidPointer(x,2); 1936 *x = snes->vec_sol_always; 1937 PetscFunctionReturn(0); 1938 } 1939 1940 #undef __FUNCT__ 1941 #define __FUNCT__ "SNESGetSolutionUpdate" 1942 /*@C 1943 SNESGetSolutionUpdate - Returns the vector where the solution update is 1944 stored. 1945 1946 Not Collective, but Vec is parallel if SNES is parallel 1947 1948 Input Parameter: 1949 . snes - the SNES context 1950 1951 Output Parameter: 1952 . x - the solution update 1953 1954 Level: advanced 1955 1956 .keywords: SNES, nonlinear, get, solution, update 1957 1958 .seealso: SNESGetSolution(), SNESGetFunction 1959 @*/ 1960 int SNESGetSolutionUpdate(SNES snes,Vec *x) 1961 { 1962 PetscFunctionBegin; 1963 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1964 PetscValidPointer(x,2); 1965 *x = snes->vec_sol_update_always; 1966 PetscFunctionReturn(0); 1967 } 1968 1969 #undef __FUNCT__ 1970 #define __FUNCT__ "SNESGetFunction" 1971 /*@C 1972 SNESGetFunction - Returns the vector where the function is stored. 1973 1974 Not Collective, but Vec is parallel if SNES is parallel 1975 1976 Input Parameter: 1977 . snes - the SNES context 1978 1979 Output Parameter: 1980 + r - the function (or PETSC_NULL) 1981 . ctx - the function context (or PETSC_NULL) 1982 - func - the function (or PETSC_NULL) 1983 1984 Level: advanced 1985 1986 .keywords: SNES, nonlinear, get, function 1987 1988 .seealso: SNESSetFunction(), SNESGetSolution() 1989 @*/ 1990 int SNESGetFunction(SNES snes,Vec *r,void **ctx,int (**func)(SNES,Vec,Vec,void*)) 1991 { 1992 PetscFunctionBegin; 1993 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 1994 if (r) *r = snes->vec_func_always; 1995 if (ctx) *ctx = snes->funP; 1996 if (func) *func = snes->computefunction; 1997 PetscFunctionReturn(0); 1998 } 1999 2000 #undef __FUNCT__ 2001 #define __FUNCT__ "SNESSetOptionsPrefix" 2002 /*@C 2003 SNESSetOptionsPrefix - Sets the prefix used for searching for all 2004 SNES options in the database. 2005 2006 Collective on SNES 2007 2008 Input Parameter: 2009 + snes - the SNES context 2010 - prefix - the prefix to prepend to all option names 2011 2012 Notes: 2013 A hyphen (-) must NOT be given at the beginning of the prefix name. 2014 The first character of all runtime options is AUTOMATICALLY the hyphen. 2015 2016 Level: advanced 2017 2018 .keywords: SNES, set, options, prefix, database 2019 2020 .seealso: SNESSetFromOptions() 2021 @*/ 2022 int SNESSetOptionsPrefix(SNES snes,const char prefix[]) 2023 { 2024 int ierr; 2025 2026 PetscFunctionBegin; 2027 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2028 ierr = PetscObjectSetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2029 ierr = KSPSetOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2030 PetscFunctionReturn(0); 2031 } 2032 2033 #undef __FUNCT__ 2034 #define __FUNCT__ "SNESAppendOptionsPrefix" 2035 /*@C 2036 SNESAppendOptionsPrefix - Appends to the prefix used for searching for all 2037 SNES options in the database. 2038 2039 Collective on SNES 2040 2041 Input Parameters: 2042 + snes - the SNES context 2043 - prefix - the prefix to prepend to all option names 2044 2045 Notes: 2046 A hyphen (-) must NOT be given at the beginning of the prefix name. 2047 The first character of all runtime options is AUTOMATICALLY the hyphen. 2048 2049 Level: advanced 2050 2051 .keywords: SNES, append, options, prefix, database 2052 2053 .seealso: SNESGetOptionsPrefix() 2054 @*/ 2055 int SNESAppendOptionsPrefix(SNES snes,const char prefix[]) 2056 { 2057 int ierr; 2058 2059 PetscFunctionBegin; 2060 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2061 ierr = PetscObjectAppendOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2062 ierr = KSPAppendOptionsPrefix(snes->ksp,prefix);CHKERRQ(ierr); 2063 PetscFunctionReturn(0); 2064 } 2065 2066 #undef __FUNCT__ 2067 #define __FUNCT__ "SNESGetOptionsPrefix" 2068 /*@C 2069 SNESGetOptionsPrefix - Sets the prefix used for searching for all 2070 SNES options in the database. 2071 2072 Not Collective 2073 2074 Input Parameter: 2075 . snes - the SNES context 2076 2077 Output Parameter: 2078 . prefix - pointer to the prefix string used 2079 2080 Notes: On the fortran side, the user should pass in a string 'prifix' of 2081 sufficient length to hold the prefix. 2082 2083 Level: advanced 2084 2085 .keywords: SNES, get, options, prefix, database 2086 2087 .seealso: SNESAppendOptionsPrefix() 2088 @*/ 2089 int SNESGetOptionsPrefix(SNES snes,char *prefix[]) 2090 { 2091 int ierr; 2092 2093 PetscFunctionBegin; 2094 PetscValidHeaderSpecific(snes,SNES_COOKIE,1); 2095 ierr = PetscObjectGetOptionsPrefix((PetscObject)snes,prefix);CHKERRQ(ierr); 2096 PetscFunctionReturn(0); 2097 } 2098 2099 2100 #undef __FUNCT__ 2101 #define __FUNCT__ "SNESRegister" 2102 /*@C 2103 SNESRegister - See SNESRegisterDynamic() 2104 2105 Level: advanced 2106 @*/ 2107 int SNESRegister(const char sname[],const char path[],const char name[],int (*function)(SNES)) 2108 { 2109 char fullname[256]; 2110 int ierr; 2111 2112 PetscFunctionBegin; 2113 ierr = PetscFListConcat(path,name,fullname);CHKERRQ(ierr); 2114 ierr = PetscFListAdd(&SNESList,sname,fullname,(void (*)(void))function);CHKERRQ(ierr); 2115 PetscFunctionReturn(0); 2116 } 2117 2118 #undef __FUNCT__ 2119 #define __FUNCT__ "SNESTestLocalMin" 2120 int SNESTestLocalMin(SNES snes) 2121 { 2122 int ierr,N,i,j; 2123 Vec u,uh,fh; 2124 PetscScalar value; 2125 PetscReal norm; 2126 2127 PetscFunctionBegin; 2128 ierr = SNESGetSolution(snes,&u);CHKERRQ(ierr); 2129 ierr = VecDuplicate(u,&uh);CHKERRQ(ierr); 2130 ierr = VecDuplicate(u,&fh);CHKERRQ(ierr); 2131 2132 /* currently only works for sequential */ 2133 ierr = PetscPrintf(PETSC_COMM_WORLD,"Testing FormFunction() for local min\n"); 2134 ierr = VecGetSize(u,&N);CHKERRQ(ierr); 2135 for (i=0; i<N; i++) { 2136 ierr = VecCopy(u,uh);CHKERRQ(ierr); 2137 ierr = PetscPrintf(PETSC_COMM_WORLD,"i = %d\n",i);CHKERRQ(ierr); 2138 for (j=-10; j<11; j++) { 2139 value = PetscSign(j)*exp(PetscAbs(j)-10.0); 2140 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2141 ierr = SNESComputeFunction(snes,uh,fh);CHKERRQ(ierr); 2142 ierr = VecNorm(fh,NORM_2,&norm);CHKERRQ(ierr); 2143 ierr = PetscPrintf(PETSC_COMM_WORLD," j norm %d %18.16e\n",j,norm);CHKERRQ(ierr); 2144 value = -value; 2145 ierr = VecSetValue(uh,i,value,ADD_VALUES);CHKERRQ(ierr); 2146 } 2147 } 2148 ierr = VecDestroy(uh);CHKERRQ(ierr); 2149 ierr = VecDestroy(fh);CHKERRQ(ierr); 2150 PetscFunctionReturn(0); 2151 } 2152