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