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