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