1 /* Defines the basic SNES object */ 2 #include <../src/snes/impls/fas/fasimpls.h> 3 4 /*MC 5 Full Approximation Scheme nonlinear multigrid solver. 6 7 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 8 9 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 10 M*/ 11 12 extern PetscErrorCode SNESDestroy_FAS(SNES snes); 13 extern PetscErrorCode SNESSetUp_FAS(SNES snes); 14 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 15 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 16 extern PetscErrorCode SNESSolve_FAS(SNES snes); 17 extern PetscErrorCode SNESReset_FAS(SNES snes); 18 19 EXTERN_C_BEGIN 20 21 #undef __FUNCT__ 22 #define __FUNCT__ "SNESCreate_FAS" 23 PetscErrorCode SNESCreate_FAS(SNES snes) 24 { 25 SNES_FAS * fas; 26 PetscErrorCode ierr; 27 28 PetscFunctionBegin; 29 snes->ops->destroy = SNESDestroy_FAS; 30 snes->ops->setup = SNESSetUp_FAS; 31 snes->ops->setfromoptions = SNESSetFromOptions_FAS; 32 snes->ops->view = SNESView_FAS; 33 snes->ops->solve = SNESSolve_FAS; 34 snes->ops->reset = SNESReset_FAS; 35 36 snes->usesksp = PETSC_FALSE; 37 snes->usespc = PETSC_FALSE; 38 39 ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 40 snes->data = (void*) fas; 41 fas->level = 0; 42 fas->levels = 1; 43 fas->n_cycles = 1; 44 fas->max_up_it = 1; 45 fas->max_down_it = 1; 46 fas->upsmooth = PETSC_NULL; 47 fas->downsmooth = PETSC_NULL; 48 fas->next = PETSC_NULL; 49 fas->interpolate = PETSC_NULL; 50 fas->restrct = PETSC_NULL; 51 fas->inject = PETSC_NULL; 52 fas->monitor = PETSC_NULL; 53 fas->usedmfornumberoflevels = PETSC_FALSE; 54 fas->useGS = PETSC_FALSE; 55 PetscFunctionReturn(0); 56 } 57 EXTERN_C_END 58 59 #undef __FUNCT__ 60 #define __FUNCT__ "SNESFASGetLevels" 61 /*@ 62 SNESFASGetLevels - Gets the number of levels to use with FAS. 63 64 Input Parameter: 65 . snes - the preconditioner context 66 67 Output parameter: 68 . levels - the number of levels 69 70 Level: advanced 71 72 .keywords: MG, get, levels, multigrid 73 74 .seealso: SNESFASSetLevels(), PCMGGetLevels() 75 @*/ 76 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 77 SNES_FAS * fas = (SNES_FAS *)snes->data; 78 PetscFunctionBegin; 79 *levels = fas->levels; 80 PetscFunctionReturn(0); 81 } 82 83 #undef __FUNCT__ 84 #define __FUNCT__ "SNESFASSetCycles" 85 /*@ 86 SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 87 complicated cycling. 88 89 Logically Collective on SNES 90 91 Input Parameters: 92 + snes - the multigrid context 93 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 94 95 Options Database Key: 96 $ -snes_fas_cycles 1 or 2 97 98 Level: advanced 99 100 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 101 102 .seealso: SNESFASSetCyclesOnLevel() 103 @*/ 104 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 105 SNES_FAS * fas = (SNES_FAS *)snes->data; 106 PetscErrorCode ierr; 107 PetscFunctionBegin; 108 fas->n_cycles = cycles; 109 if (fas->next) { 110 ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 111 } 112 PetscFunctionReturn(0); 113 } 114 115 #undef __FUNCT__ 116 #define __FUNCT__ "SNESFASSetCyclesOnLevel" 117 /*@ 118 SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 119 complicated cycling. 120 121 Logically Collective on SNES 122 123 Input Parameters: 124 + snes - the multigrid context 125 . level - the level to set the number of cycles on 126 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 127 128 Level: advanced 129 130 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 131 132 .seealso: SNESFASSetCycles() 133 @*/ 134 PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 135 SNES_FAS * fas = (SNES_FAS *)snes->data; 136 PetscInt top_level = fas->level,i; 137 138 PetscFunctionBegin; 139 if (level > top_level) 140 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 141 /* get to the correct level */ 142 for (i = fas->level; i > level; i--) { 143 fas = (SNES_FAS *)fas->next->data; 144 } 145 if (fas->level != level) 146 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 147 fas->n_cycles = cycles; 148 PetscFunctionReturn(0); 149 } 150 151 #undef __FUNCT__ 152 #define __FUNCT__ "SNESFASSetGS" 153 /*@ 154 SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 155 Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 156 and nonlinear preconditioners. 157 158 Logically Collective on SNES 159 160 Input Parameters: 161 + snes - the multigrid context 162 . gsfunc - the nonlinear smoother function 163 . ctx - the user context for the nonlinear smoother 164 - use_gs - whether to use the nonlinear smoother or not 165 166 Level: advanced 167 168 .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 169 170 .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 171 @*/ 172 PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 173 PetscErrorCode ierr = 0; 174 SNES_FAS *fas = (SNES_FAS *)snes->data; 175 PetscFunctionBegin; 176 177 /* use or don't use it according to user wishes*/ 178 fas->useGS = use_gs; 179 if (gsfunc) { 180 ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 181 /* push the provided GS up the tree */ 182 if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 183 } else if (snes->ops->computegs) { 184 /* assume that the user has set the GS solver at this level */ 185 if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 186 } else if (use_gs) { 187 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 188 } 189 PetscFunctionReturn(0); 190 } 191 192 #undef __FUNCT__ 193 #define __FUNCT__ "SNESFASSetGSOnLevel" 194 /*@ 195 SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 196 197 Logically Collective on SNES 198 199 Input Parameters: 200 + snes - the multigrid context 201 . level - the level to set the nonlinear smoother on 202 . gsfunc - the nonlinear smoother function 203 . ctx - the user context for the nonlinear smoother 204 - use_gs - whether to use the nonlinear smoother or not 205 206 Level: advanced 207 208 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 209 210 .seealso: SNESSetGS(), SNESFASSetGS() 211 @*/ 212 PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 213 SNES_FAS *fas = (SNES_FAS *)snes->data; 214 PetscErrorCode ierr; 215 PetscInt top_level = fas->level,i; 216 SNES cur_snes = snes; 217 PetscFunctionBegin; 218 if (level > top_level) 219 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 220 /* get to the correct level */ 221 for (i = fas->level; i > level; i--) { 222 fas = (SNES_FAS *)fas->next->data; 223 cur_snes = fas->next; 224 } 225 if (fas->level != level) 226 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 227 fas->useGS = use_gs; 228 if (gsfunc) { 229 ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 230 } 231 PetscFunctionReturn(0); 232 } 233 234 #undef __FUNCT__ 235 #define __FUNCT__ "SNESFASGetSNES" 236 /*@ 237 SNESFASGetSNES - Gets the SNES corresponding to a particular 238 level of the FAS hierarchy. 239 240 Input Parameters: 241 + snes - the multigrid context 242 level - the level to get 243 - lsnes - whether to use the nonlinear smoother or not 244 245 Level: advanced 246 247 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 248 249 .seealso: SNESFASSetLevels(), SNESFASGetLevels() 250 @*/ 251 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 252 SNES_FAS * fas = (SNES_FAS *)snes->data; 253 PetscInt levels = fas->level; 254 PetscInt i; 255 PetscFunctionBegin; 256 *lsnes = snes; 257 if (fas->level < level) { 258 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 259 } 260 if (level > levels - 1) { 261 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 262 } 263 for (i = fas->level; i > level; i--) { 264 *lsnes = fas->next; 265 fas = (SNES_FAS *)(*lsnes)->data; 266 } 267 if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 268 PetscFunctionReturn(0); 269 } 270 271 #undef __FUNCT__ 272 #define __FUNCT__ "SNESFASSetLevels" 273 /*@ 274 SNESFASSetLevels - Sets the number of levels to use with FAS. 275 Must be called before any other FAS routine. 276 277 Input Parameters: 278 + snes - the snes context 279 . levels - the number of levels 280 - comms - optional communicators for each level; this is to allow solving the coarser 281 problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 282 Fortran. 283 284 Level: intermediate 285 286 Notes: 287 If the number of levels is one then the multigrid uses the -fas_levels prefix 288 for setting the level options rather than the -fas_coarse prefix. 289 290 .keywords: FAS, MG, set, levels, multigrid 291 292 .seealso: SNESFASGetLevels() 293 @*/ 294 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 295 PetscErrorCode ierr; 296 PetscInt i; 297 SNES_FAS * fas = (SNES_FAS *)snes->data; 298 MPI_Comm comm; 299 PetscFunctionBegin; 300 comm = ((PetscObject)snes)->comm; 301 if (levels == fas->levels) { 302 if (!comms) 303 PetscFunctionReturn(0); 304 } 305 /* user has changed the number of levels; reset */ 306 ierr = SNESReset(snes);CHKERRQ(ierr); 307 /* destroy any coarser levels if necessary */ 308 if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 309 fas->next = PETSC_NULL; 310 /* setup the finest level */ 311 for (i = levels - 1; i >= 0; i--) { 312 if (comms) comm = comms[i]; 313 fas->level = i; 314 fas->levels = levels; 315 fas->next = PETSC_NULL; 316 if (i > 0) { 317 ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 318 ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 319 ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 320 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 321 fas = (SNES_FAS *)fas->next->data; 322 } 323 } 324 /* by default set the GS smoothers up the chain if one exists on the finest level */ 325 if (snes->ops->computegs) 326 ierr = SNESFASSetGS(snes, snes->ops->computegs, snes->gsP, fas->useGS);CHKERRQ(ierr); 327 PetscFunctionReturn(0); 328 } 329 330 #undef __FUNCT__ 331 #define __FUNCT__ "SNESFASSetNumberSmoothUp" 332 /*@ 333 SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 334 use on all levels. 335 336 Logically Collective on PC 337 338 Input Parameters: 339 + snes - the multigrid context 340 - n - the number of smoothing steps 341 342 Options Database Key: 343 . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 344 345 Level: advanced 346 347 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 348 349 .seealso: SNESFASSetNumberSmoothDown() 350 @*/ 351 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 352 SNES_FAS * fas = (SNES_FAS *)snes->data; 353 PetscErrorCode ierr = 0; 354 PetscFunctionBegin; 355 fas->max_up_it = n; 356 if (fas->next) { 357 ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 358 } 359 PetscFunctionReturn(0); 360 } 361 362 #undef __FUNCT__ 363 #define __FUNCT__ "SNESFASSetNumberSmoothDown" 364 /*@ 365 SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 366 use on all levels. 367 368 Logically Collective on PC 369 370 Input Parameters: 371 + snes - the multigrid context 372 - n - the number of smoothing steps 373 374 Options Database Key: 375 . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 376 377 Level: advanced 378 379 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 380 381 .seealso: SNESFASSetNumberSmoothUp() 382 @*/ 383 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 384 SNES_FAS * fas = (SNES_FAS *)snes->data; 385 PetscErrorCode ierr = 0; 386 PetscFunctionBegin; 387 fas->max_down_it = n; 388 if (fas->next) { 389 ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 390 } 391 PetscFunctionReturn(0); 392 } 393 394 #undef __FUNCT__ 395 #define __FUNCT__ "SNESFASSetInterpolation" 396 /*@ 397 SNESFASSetInterpolation - Sets the function to be used to calculate the 398 interpolation from l-1 to the lth level 399 400 Logically Collective on PC and Mat 401 402 Input Parameters: 403 + snes - the multigrid context 404 . mat - the interpolation operator 405 - level - the level (0 is coarsest) to supply [do not supply 0] 406 407 Level: advanced 408 409 Notes: 410 Usually this is the same matrix used also to set the restriction 411 for the same level. 412 413 One can pass in the interpolation matrix or its transpose; PETSc figures 414 out from the matrix size which one it is. 415 416 .keywords: FAS, multigrid, set, interpolate, level 417 418 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale() 419 @*/ 420 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 421 SNES_FAS * fas = (SNES_FAS *)snes->data; 422 PetscInt top_level = fas->level,i; 423 424 PetscFunctionBegin; 425 if (level > top_level) 426 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 427 /* get to the correct level */ 428 for (i = fas->level; i > level; i--) { 429 fas = (SNES_FAS *)fas->next->data; 430 } 431 if (fas->level != level) 432 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 433 fas->interpolate = mat; 434 PetscFunctionReturn(0); 435 } 436 437 #undef __FUNCT__ 438 #define __FUNCT__ "SNESFASSetRestriction" 439 /*@ 440 SNESFASSetRestriction - Sets the function to be used to restrict the defect 441 from level l to l-1. 442 443 Logically Collective on SNES and Mat 444 445 Input Parameters: 446 + snes - the multigrid context 447 . mat - the restriction matrix 448 - level - the level (0 is coarsest) to supply [Do not supply 0] 449 450 Level: advanced 451 452 Notes: 453 Usually this is the same matrix used also to set the interpolation 454 for the same level. 455 456 One can pass in the interpolation matrix or its transpose; PETSc figures 457 out from the matrix size which one it is. 458 459 If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 460 is used. 461 462 .keywords: FAS, MG, set, multigrid, restriction, level 463 464 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 465 @*/ 466 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 467 SNES_FAS * fas = (SNES_FAS *)snes->data; 468 PetscInt top_level = fas->level,i; 469 470 PetscFunctionBegin; 471 if (level > top_level) 472 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 473 /* get to the correct level */ 474 for (i = fas->level; i > level; i--) { 475 fas = (SNES_FAS *)fas->next->data; 476 } 477 if (fas->level != level) 478 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 479 fas->restrct = mat; 480 PetscFunctionReturn(0); 481 } 482 483 484 #undef __FUNCT__ 485 #define __FUNCT__ "SNESFASSetRScale" 486 /*@ 487 SNESFASSetRScale - Sets the scaling factor of the restriction 488 operator from level l to l-1. 489 490 Logically Collective on SNES and Mat 491 492 Input Parameters: 493 + snes - the multigrid context 494 . rscale - the restriction scaling 495 - level - the level (0 is coarsest) to supply [Do not supply 0] 496 497 Level: advanced 498 499 Notes: 500 This is only used in the case that the injection is not set. 501 502 .keywords: FAS, MG, set, multigrid, restriction, level 503 504 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 505 @*/ 506 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 507 SNES_FAS * fas = (SNES_FAS *)snes->data; 508 PetscInt top_level = fas->level,i; 509 510 PetscFunctionBegin; 511 if (level > top_level) 512 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 513 /* get to the correct level */ 514 for (i = fas->level; i > level; i--) { 515 fas = (SNES_FAS *)fas->next->data; 516 } 517 if (fas->level != level) 518 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 519 fas->rscale = rscale; 520 PetscFunctionReturn(0); 521 } 522 523 524 #undef __FUNCT__ 525 #define __FUNCT__ "SNESFASSetInjection" 526 /*@ 527 SNESFASSetInjection - Sets the function to be used to inject the solution 528 from level l to l-1. 529 530 Logically Collective on SNES and Mat 531 532 Input Parameters: 533 + snes - the multigrid context 534 . mat - the restriction matrix 535 - level - the level (0 is coarsest) to supply [Do not supply 0] 536 537 Level: advanced 538 539 Notes: 540 If you do not set this, the restriction and rscale is used to 541 project the solution instead. 542 543 .keywords: FAS, MG, set, multigrid, restriction, level 544 545 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 546 @*/ 547 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 548 SNES_FAS * fas = (SNES_FAS *)snes->data; 549 PetscInt top_level = fas->level,i; 550 551 PetscFunctionBegin; 552 if (level > top_level) 553 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 554 /* get to the correct level */ 555 for (i = fas->level; i > level; i--) { 556 fas = (SNES_FAS *)fas->next->data; 557 } 558 if (fas->level != level) 559 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 560 fas->inject = mat; 561 PetscFunctionReturn(0); 562 } 563 564 #undef __FUNCT__ 565 #define __FUNCT__ "SNESReset_FAS" 566 PetscErrorCode SNESReset_FAS(SNES snes) 567 { 568 PetscErrorCode ierr = 0; 569 SNES_FAS * fas = (SNES_FAS *)snes->data; 570 571 PetscFunctionBegin; 572 /* destroy local data created in SNESSetup_FAS */ 573 #if 0 574 /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */ 575 #endif 576 if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr); 577 #if 0 578 #endif 579 PetscFunctionReturn(0); 580 } 581 582 #undef __FUNCT__ 583 #define __FUNCT__ "SNESDestroy_FAS" 584 PetscErrorCode SNESDestroy_FAS(SNES snes) 585 { 586 SNES_FAS * fas = (SNES_FAS *)snes->data; 587 PetscErrorCode ierr; 588 589 PetscFunctionBegin; 590 /* recursively resets and then destroys */ 591 ierr = SNESReset_FAS(snes);CHKERRQ(ierr); 592 if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 593 if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 594 if (fas->inject) { 595 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 596 } 597 if (fas->interpolate == fas->restrct) { 598 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 599 fas->restrct = PETSC_NULL; 600 } else { 601 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 602 if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 603 } 604 if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 605 if (snes->work) ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr); 606 if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 607 ierr = PetscFree(fas);CHKERRQ(ierr); 608 609 PetscFunctionReturn(0); 610 } 611 612 #undef __FUNCT__ 613 #define __FUNCT__ "SNESSetUp_FAS" 614 PetscErrorCode SNESSetUp_FAS(SNES snes) 615 { 616 SNES_FAS *fas = (SNES_FAS *) snes->data,*tmp; 617 PetscErrorCode ierr; 618 VecScatter injscatter; 619 PetscInt dm_levels; 620 621 622 PetscFunctionBegin; 623 624 if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 625 ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 626 dm_levels++; 627 if (dm_levels > fas->levels) { 628 ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 629 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 630 } 631 } 632 633 if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */ 634 635 /* should call the SNESSetFromOptions() only when appropriate */ 636 tmp = fas; 637 while (tmp) { 638 if (tmp->upsmooth) {ierr = SNESSetFromOptions(tmp->upsmooth);CHKERRQ(ierr);} 639 if (tmp->downsmooth) {ierr = SNESSetFromOptions(tmp->downsmooth);CHKERRQ(ierr);} 640 tmp = tmp->next ? (SNES_FAS*) tmp->next->data : 0; 641 } 642 643 /* gets the solver ready for solution */ 644 if (fas->next) { 645 if (snes->ops->computegs) { 646 ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 647 } 648 } 649 if (snes->dm) { 650 /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 651 if (fas->next) { 652 /* for now -- assume the DM and the evaluation functions have been set externally */ 653 if (!fas->next->dm) { 654 ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 655 ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 656 } 657 /* set the interpolation and restriction from the DM */ 658 if (!fas->interpolate) { 659 ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 660 fas->restrct = fas->interpolate; 661 } 662 /* set the injection from the DM */ 663 if (!fas->inject) { 664 ierr = DMGetInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 665 ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 666 ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 667 } 668 } 669 /* set the DMs of the pre and post-smoothers here */ 670 if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 671 if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 672 } 673 674 if (fas->next) { 675 /* gotta set up the solution vector for this to work */ 676 if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 677 ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 678 } 679 /* got to set them all up at once */ 680 PetscFunctionReturn(0); 681 } 682 683 #undef __FUNCT__ 684 #define __FUNCT__ "SNESSetFromOptions_FAS" 685 PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 686 { 687 SNES_FAS *fas = (SNES_FAS *) snes->data; 688 PetscInt levels = 1; 689 PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg; 690 PetscErrorCode ierr; 691 const char * def_smooth = SNESNRICHARDSON; 692 char pre_type[256]; 693 char post_type[256]; 694 char monfilename[PETSC_MAX_PATH_LEN]; 695 696 PetscFunctionBegin; 697 ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 698 699 /* number of levels -- only process on the finest level */ 700 if (fas->levels - 1 == fas->level) { 701 ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 702 if (!flg && snes->dm) { 703 ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 704 levels++; 705 fas->usedmfornumberoflevels = PETSC_TRUE; 706 } 707 ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 708 } 709 710 /* type of pre and/or post smoothers -- set both at once */ 711 ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 712 ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 713 ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 714 if (smoothflg) { 715 ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 716 } else { 717 ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr); 718 ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr); 719 } 720 721 /* options for the number of preconditioning cycles and cycle type */ 722 ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 723 ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 724 ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 725 726 ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 727 728 ierr = PetscOptionsBool("-snes_fas_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr); 729 730 /* other options for the coarsest level */ 731 if (fas->level == 0) { 732 ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 733 ierr = PetscOptionsBool("-snes_fas_coarse_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr); 734 } 735 736 ierr = PetscOptionsTail();CHKERRQ(ierr); 737 /* setup from the determined types if there is no pointwise procedure or smoother defined */ 738 739 if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !fas->useGS)) { 740 const char *prefix; 741 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 742 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 743 ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 744 if (fas->level || (fas->levels == 1)) { 745 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr); 746 } else { 747 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 748 } 749 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 750 ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 751 if (snes->ops->computefunction) { 752 ierr = SNESSetFunction(fas->downsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 753 } 754 } 755 756 if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !fas->useGS)) { 757 const char *prefix; 758 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 759 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 760 ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 761 ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr); 762 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 763 ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 764 if (snes->ops->computefunction) { 765 ierr = SNESSetFunction(fas->upsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr); 766 } 767 } 768 if (fas->upsmooth) { 769 ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 770 } 771 772 if (fas->downsmooth) { 773 ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 774 } 775 776 if (monflg) { 777 fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 778 /* set the monitors for the upsmoother and downsmoother */ 779 if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 780 if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 781 } 782 783 /* recursive option setting for the smoothers */ 784 if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);} 785 PetscFunctionReturn(0); 786 } 787 788 #undef __FUNCT__ 789 #define __FUNCT__ "SNESView_FAS" 790 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 791 { 792 SNES_FAS *fas = (SNES_FAS *) snes->data; 793 PetscBool iascii; 794 PetscErrorCode ierr; 795 PetscInt levels = fas->levels; 796 PetscInt i; 797 798 PetscFunctionBegin; 799 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 800 if (iascii) { 801 ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 802 ierr = PetscViewerASCIIPushTab(viewer); 803 for (i = levels - 1; i >= 0; i--) { 804 ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 805 if (fas->upsmooth) { 806 ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 807 ierr = PetscViewerASCIIPushTab(viewer); 808 ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 809 ierr = PetscViewerASCIIPopTab(viewer); 810 } else { 811 ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 812 } 813 if (fas->downsmooth) { 814 ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 815 ierr = PetscViewerASCIIPushTab(viewer); 816 ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 817 ierr = PetscViewerASCIIPopTab(viewer); 818 } else { 819 ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 820 } 821 if (fas->useGS) { 822 ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n", fas->level);CHKERRQ(ierr); 823 } 824 if (fas->next) fas = (SNES_FAS *)fas->next->data; 825 } 826 ierr = PetscViewerASCIIPopTab(viewer); 827 } else { 828 SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 829 } 830 PetscFunctionReturn(0); 831 } 832 833 /* 834 835 Defines the FAS cycle as: 836 837 fine problem: F(x) = 0 838 coarse problem: F^c(x) = b^c 839 840 b^c = F^c(I^c_fx^f - I^c_fF(x)) 841 842 correction: 843 844 x = x + I(x^c - Rx) 845 846 */ 847 848 #undef __FUNCT__ 849 #define __FUNCT__ "FASCycle_Private" 850 PetscErrorCode FASCycle_Private(SNES snes, Vec X) { 851 852 PetscErrorCode ierr; 853 Vec X_c, Xo_c, F_c, B_c,F,B; 854 SNES_FAS * fas = (SNES_FAS *)snes->data; 855 PetscInt i, k; 856 PetscReal fnorm; 857 858 PetscFunctionBegin; 859 F = snes->vec_func; 860 B = snes->vec_rhs; 861 /* pre-smooth -- just update using the pre-smoother */ 862 if (fas->downsmooth) { 863 ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 864 } else if (snes->ops->computegs) { 865 if (fas->monitor) { 866 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 867 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 868 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 869 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 870 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 871 } 872 for (k = 0; k < fas->max_up_it; k++) { 873 ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 874 if (fas->monitor) { 875 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 876 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 877 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 878 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 879 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 880 } 881 } 882 } else if (snes->pc) { 883 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 884 } 885 if (fas->next) { 886 for (i = 0; i < fas->n_cycles; i++) { 887 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 888 889 X_c = fas->next->vec_sol; 890 Xo_c = fas->next->work[0]; 891 F_c = fas->next->vec_func; 892 B_c = fas->next->work[1]; 893 894 /* inject the solution */ 895 if (fas->inject) { 896 ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 897 } else { 898 ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 899 ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 900 } 901 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 902 903 /* restrict the defect */ 904 ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 905 906 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 907 fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 908 ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 909 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 910 911 /* set initial guess of the coarse problem to the projected fine solution */ 912 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 913 914 /* recurse to the next level */ 915 fas->next->vec_rhs = B_c; 916 ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); 917 fas->next->vec_rhs = PETSC_NULL; 918 919 /* correct as x <- x + I(x^c - Rx)*/ 920 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 921 ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 922 } 923 } 924 /* up-smooth -- just update using the down-smoother */ 925 if (fas->level != 0) { 926 if (fas->upsmooth) { 927 ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 928 } else if (snes->ops->computegs) { 929 if (fas->monitor) { 930 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 931 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 932 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 933 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 934 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 935 } 936 for (k = 0; k < fas->max_down_it; k++) { 937 ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 938 if (fas->monitor) { 939 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 940 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 941 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 942 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 943 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 944 } 945 } 946 } else if (snes->pc) { 947 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 948 } 949 } 950 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 951 952 PetscFunctionReturn(0); 953 } 954 955 #undef __FUNCT__ 956 #define __FUNCT__ "SNESSolve_FAS" 957 958 PetscErrorCode SNESSolve_FAS(SNES snes) 959 { 960 PetscErrorCode ierr; 961 PetscInt i, maxits; 962 Vec X, B,F; 963 PetscReal fnorm; 964 PetscFunctionBegin; 965 maxits = snes->max_its; /* maximum number of iterations */ 966 snes->reason = SNES_CONVERGED_ITERATING; 967 X = snes->vec_sol; 968 B = snes->vec_rhs; 969 F = snes->vec_func; 970 971 /*norm setup */ 972 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 973 snes->iter = 0; 974 snes->norm = 0.; 975 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 976 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 977 if (snes->domainerror) { 978 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 979 PetscFunctionReturn(0); 980 } 981 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 982 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 983 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 984 snes->norm = fnorm; 985 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 986 SNESLogConvHistory(snes,fnorm,0); 987 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 988 989 /* set parameter for default relative tolerance convergence test */ 990 snes->ttol = fnorm*snes->rtol; 991 /* test convergence */ 992 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 993 if (snes->reason) PetscFunctionReturn(0); 994 for (i = 0; i < maxits; i++) { 995 /* Call general purpose update function */ 996 997 if (snes->ops->update) { 998 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 999 } 1000 ierr = FASCycle_Private(snes, X);CHKERRQ(ierr); 1001 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1002 /* Monitor convergence */ 1003 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1004 snes->iter = i+1; 1005 snes->norm = fnorm; 1006 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1007 SNESLogConvHistory(snes,snes->norm,0); 1008 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1009 /* Test for convergence */ 1010 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1011 if (snes->reason) break; 1012 } 1013 if (i == maxits) { 1014 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1015 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1016 } 1017 PetscFunctionReturn(0); 1018 } 1019