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 extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 19 20 EXTERN_C_BEGIN 21 22 #undef __FUNCT__ 23 #define __FUNCT__ "SNESCreate_FAS" 24 PetscErrorCode SNESCreate_FAS(SNES snes) 25 { 26 SNES_FAS * fas; 27 PetscErrorCode ierr; 28 29 PetscFunctionBegin; 30 snes->ops->destroy = SNESDestroy_FAS; 31 snes->ops->setup = SNESSetUp_FAS; 32 snes->ops->setfromoptions = SNESSetFromOptions_FAS; 33 snes->ops->view = SNESView_FAS; 34 snes->ops->solve = SNESSolve_FAS; 35 snes->ops->reset = SNESReset_FAS; 36 37 snes->usesksp = PETSC_FALSE; 38 snes->usespc = PETSC_FALSE; 39 40 ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 41 snes->data = (void*) fas; 42 fas->level = 0; 43 fas->levels = 1; 44 fas->n_cycles = 1; 45 fas->max_up_it = 1; 46 fas->max_down_it = 1; 47 fas->upsmooth = PETSC_NULL; 48 fas->downsmooth = PETSC_NULL; 49 fas->next = PETSC_NULL; 50 fas->previous = PETSC_NULL; 51 fas->interpolate = PETSC_NULL; 52 fas->restrct = PETSC_NULL; 53 fas->inject = PETSC_NULL; 54 fas->monitor = PETSC_NULL; 55 fas->usedmfornumberoflevels = PETSC_FALSE; 56 PetscFunctionReturn(0); 57 } 58 EXTERN_C_END 59 60 #undef __FUNCT__ 61 #define __FUNCT__ "SNESFASGetLevels" 62 /*@ 63 SNESFASGetLevels - Gets the number of levels in a FAS. 64 65 Input Parameter: 66 . snes - the preconditioner context 67 68 Output parameter: 69 . levels - the number of levels 70 71 Level: advanced 72 73 .keywords: MG, get, levels, multigrid 74 75 .seealso: SNESFASSetLevels(), PCMGGetLevels() 76 @*/ 77 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 78 SNES_FAS * fas = (SNES_FAS *)snes->data; 79 PetscFunctionBegin; 80 *levels = fas->levels; 81 PetscFunctionReturn(0); 82 } 83 84 #undef __FUNCT__ 85 #define __FUNCT__ "SNESFASSetCycles" 86 /*@ 87 SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 88 complicated cycling. 89 90 Logically Collective on SNES 91 92 Input Parameters: 93 + snes - the multigrid context 94 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 95 96 Options Database Key: 97 $ -snes_fas_cycles 1 or 2 98 99 Level: advanced 100 101 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 102 103 .seealso: SNESFASSetCyclesOnLevel() 104 @*/ 105 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 106 SNES_FAS * fas = (SNES_FAS *)snes->data; 107 PetscErrorCode ierr; 108 PetscFunctionBegin; 109 fas->n_cycles = cycles; 110 if (fas->next) { 111 ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 112 } 113 PetscFunctionReturn(0); 114 } 115 116 #undef __FUNCT__ 117 #define __FUNCT__ "SNESFASSetCyclesOnLevel" 118 /*@ 119 SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 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 /*@C 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 snes->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 /*@C 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 snes->usegs = use_gs; 228 if (gsfunc) { 229 ierr = SNESSetGS(cur_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 /*@C 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 SNES prevsnes; 299 MPI_Comm comm; 300 PetscFunctionBegin; 301 comm = ((PetscObject)snes)->comm; 302 if (levels == fas->levels) { 303 if (!comms) 304 PetscFunctionReturn(0); 305 } 306 /* user has changed the number of levels; reset */ 307 ierr = SNESReset(snes);CHKERRQ(ierr); 308 /* destroy any coarser levels if necessary */ 309 if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 310 fas->next = PETSC_NULL; 311 fas->previous = PETSC_NULL; 312 prevsnes = snes; 313 /* setup the finest level */ 314 for (i = levels - 1; i >= 0; i--) { 315 if (comms) comm = comms[i]; 316 fas->level = i; 317 fas->levels = levels; 318 fas->next = PETSC_NULL; 319 if (i > 0) { 320 ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 321 ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 322 ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 323 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 324 ((SNES_FAS *)fas->next->data)->previous = prevsnes; 325 prevsnes = fas->next; 326 fas = (SNES_FAS *)prevsnes->data; 327 } 328 } 329 PetscFunctionReturn(0); 330 } 331 332 #undef __FUNCT__ 333 #define __FUNCT__ "SNESFASSetNumberSmoothUp" 334 /*@ 335 SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 336 use on all levels. 337 338 Logically Collective on PC 339 340 Input Parameters: 341 + snes - the multigrid context 342 - n - the number of smoothing steps 343 344 Options Database Key: 345 . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 346 347 Level: advanced 348 349 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 350 351 .seealso: SNESFASSetNumberSmoothDown() 352 @*/ 353 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 354 SNES_FAS * fas = (SNES_FAS *)snes->data; 355 PetscErrorCode ierr = 0; 356 PetscFunctionBegin; 357 fas->max_up_it = n; 358 if (fas->next) { 359 ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 360 } 361 PetscFunctionReturn(0); 362 } 363 364 #undef __FUNCT__ 365 #define __FUNCT__ "SNESFASSetNumberSmoothDown" 366 /*@ 367 SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 368 use on all levels. 369 370 Logically Collective on PC 371 372 Input Parameters: 373 + snes - the multigrid context 374 - n - the number of smoothing steps 375 376 Options Database Key: 377 . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 378 379 Level: advanced 380 381 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 382 383 .seealso: SNESFASSetNumberSmoothUp() 384 @*/ 385 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 386 SNES_FAS * fas = (SNES_FAS *)snes->data; 387 PetscErrorCode ierr = 0; 388 PetscFunctionBegin; 389 fas->max_down_it = n; 390 if (fas->next) { 391 ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 392 } 393 PetscFunctionReturn(0); 394 } 395 396 #undef __FUNCT__ 397 #define __FUNCT__ "SNESFASSetInterpolation" 398 /*@ 399 SNESFASSetInterpolation - Sets the function to be used to calculate the 400 interpolation from l-1 to the lth level 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 Input Parameters: 444 + snes - the multigrid context 445 . mat - the restriction matrix 446 - level - the level (0 is coarsest) to supply [Do not supply 0] 447 448 Level: advanced 449 450 Notes: 451 Usually this is the same matrix used also to set the interpolation 452 for the same level. 453 454 One can pass in the interpolation matrix or its transpose; PETSc figures 455 out from the matrix size which one it is. 456 457 If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 458 is used. 459 460 .keywords: FAS, MG, set, multigrid, restriction, level 461 462 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 463 @*/ 464 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 465 SNES_FAS * fas = (SNES_FAS *)snes->data; 466 PetscInt top_level = fas->level,i; 467 468 PetscFunctionBegin; 469 if (level > top_level) 470 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 471 /* get to the correct level */ 472 for (i = fas->level; i > level; i--) { 473 fas = (SNES_FAS *)fas->next->data; 474 } 475 if (fas->level != level) 476 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 477 fas->restrct = mat; 478 PetscFunctionReturn(0); 479 } 480 481 482 #undef __FUNCT__ 483 #define __FUNCT__ "SNESFASSetRScale" 484 /*@ 485 SNESFASSetRScale - Sets the scaling factor of the restriction 486 operator from level l to l-1. 487 488 Input Parameters: 489 + snes - the multigrid context 490 . rscale - the restriction scaling 491 - level - the level (0 is coarsest) to supply [Do not supply 0] 492 493 Level: advanced 494 495 Notes: 496 This is only used in the case that the injection is not set. 497 498 .keywords: FAS, MG, set, multigrid, restriction, level 499 500 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 501 @*/ 502 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 503 SNES_FAS * fas = (SNES_FAS *)snes->data; 504 PetscInt top_level = fas->level,i; 505 506 PetscFunctionBegin; 507 if (level > top_level) 508 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 509 /* get to the correct level */ 510 for (i = fas->level; i > level; i--) { 511 fas = (SNES_FAS *)fas->next->data; 512 } 513 if (fas->level != level) 514 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 515 fas->rscale = rscale; 516 PetscFunctionReturn(0); 517 } 518 519 520 #undef __FUNCT__ 521 #define __FUNCT__ "SNESFASSetInjection" 522 /*@ 523 SNESFASSetInjection - Sets the function to be used to inject the solution 524 from level l to l-1. 525 526 Input Parameters: 527 + snes - the multigrid context 528 . mat - the restriction matrix 529 - level - the level (0 is coarsest) to supply [Do not supply 0] 530 531 Level: advanced 532 533 Notes: 534 If you do not set this, the restriction and rscale is used to 535 project the solution instead. 536 537 .keywords: FAS, MG, set, multigrid, restriction, level 538 539 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 540 @*/ 541 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 542 SNES_FAS * fas = (SNES_FAS *)snes->data; 543 PetscInt top_level = fas->level,i; 544 545 PetscFunctionBegin; 546 if (level > top_level) 547 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 548 /* get to the correct level */ 549 for (i = fas->level; i > level; i--) { 550 fas = (SNES_FAS *)fas->next->data; 551 } 552 if (fas->level != level) 553 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 554 fas->inject = mat; 555 PetscFunctionReturn(0); 556 } 557 558 #undef __FUNCT__ 559 #define __FUNCT__ "SNESReset_FAS" 560 PetscErrorCode SNESReset_FAS(SNES snes) 561 { 562 PetscErrorCode ierr = 0; 563 SNES_FAS * fas = (SNES_FAS *)snes->data; 564 565 PetscFunctionBegin; 566 if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 567 if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 568 if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 569 PetscFunctionReturn(0); 570 } 571 572 #undef __FUNCT__ 573 #define __FUNCT__ "SNESDestroy_FAS" 574 PetscErrorCode SNESDestroy_FAS(SNES snes) 575 { 576 SNES_FAS * fas = (SNES_FAS *)snes->data; 577 PetscErrorCode ierr = 0; 578 579 PetscFunctionBegin; 580 /* recursively resets and then destroys */ 581 if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 582 if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 583 if (fas->inject) { 584 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 585 } 586 if (fas->interpolate == fas->restrct) { 587 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 588 fas->restrct = PETSC_NULL; 589 } else { 590 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 591 if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 592 } 593 if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 594 if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 595 ierr = PetscFree(fas);CHKERRQ(ierr); 596 597 PetscFunctionReturn(0); 598 } 599 600 #undef __FUNCT__ 601 #define __FUNCT__ "SNESSetUp_FAS" 602 PetscErrorCode SNESSetUp_FAS(SNES snes) 603 { 604 SNES_FAS *fas = (SNES_FAS *) snes->data; 605 PetscErrorCode ierr; 606 VecScatter injscatter; 607 PetscInt dm_levels; 608 609 PetscFunctionBegin; 610 611 if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 612 ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 613 dm_levels++; 614 if (dm_levels > fas->levels) { 615 ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 616 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 617 } 618 } 619 620 ierr = SNESDefaultGetWork(snes, 1);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 621 622 /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 623 if (fas->upsmooth) { 624 ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 625 if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 626 ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 627 } 628 if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 629 ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 630 } 631 if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 632 ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 633 } 634 } 635 if (fas->downsmooth) { 636 ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 637 if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 638 ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 639 } 640 if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 641 ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 642 } 643 if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 644 ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 645 } 646 } 647 /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 648 if (fas->next) { 649 if (fas->galerkin) { 650 ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 651 } else { 652 if (snes->ops->computefunction && !fas->next->ops->computefunction) { 653 ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 654 } 655 if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 656 ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 657 } 658 if (snes->ops->computegs && !fas->next->ops->computegs) { 659 ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 660 } 661 } 662 } 663 if (snes->dm) { 664 /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 665 if (fas->next) { 666 /* for now -- assume the DM and the evaluation functions have been set externally */ 667 if (!fas->next->dm) { 668 ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 669 ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 670 } 671 /* set the interpolation and restriction from the DM */ 672 if (!fas->interpolate) { 673 ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 674 fas->restrct = fas->interpolate; 675 } 676 /* set the injection from the DM */ 677 if (!fas->inject) { 678 ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 679 ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 680 ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 681 } 682 } 683 /* set the DMs of the pre and post-smoothers here */ 684 if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 685 if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 686 } 687 688 /* setup FAS work vectors */ 689 if (fas->galerkin) { 690 ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 691 ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 692 } 693 694 if (fas->next) { 695 /* gotta set up the solution vector for this to work */ 696 if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 697 if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);} 698 ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 699 } 700 /* got to set them all up at once */ 701 PetscFunctionReturn(0); 702 } 703 704 #undef __FUNCT__ 705 #define __FUNCT__ "SNESSetFromOptions_FAS" 706 PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 707 { 708 SNES_FAS *fas = (SNES_FAS *) snes->data; 709 PetscInt levels = 1; 710 PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg; 711 PetscErrorCode ierr; 712 const char *def_smooth = SNESNRICHARDSON; 713 char pre_type[256]; 714 char post_type[256]; 715 char monfilename[PETSC_MAX_PATH_LEN]; 716 717 PetscFunctionBegin; 718 ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 719 720 /* number of levels -- only process on the finest level */ 721 if (fas->levels - 1 == fas->level) { 722 ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 723 if (!flg && snes->dm) { 724 ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 725 levels++; 726 fas->usedmfornumberoflevels = PETSC_TRUE; 727 } 728 ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 729 } 730 731 /* type of pre and/or post smoothers -- set both at once */ 732 ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 733 ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 734 ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 735 if (smoothflg) { 736 ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 737 } else { 738 ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr); 739 ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr); 740 } 741 742 /* options for the number of preconditioning cycles and cycle type */ 743 ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 744 ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 745 ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 746 747 ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 748 749 ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 750 751 /* other options for the coarsest level */ 752 if (fas->level == 0) { 753 ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 754 } 755 756 ierr = PetscOptionsTail();CHKERRQ(ierr); 757 /* setup from the determined types if there is no pointwise procedure or smoother defined */ 758 759 if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->usegs)) { 760 const char *prefix; 761 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 762 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 763 ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 764 if (fas->level || (fas->levels == 1)) { 765 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr); 766 } else { 767 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 768 } 769 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 770 ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 771 } 772 773 if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->usegs)) { 774 const char *prefix; 775 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 776 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 777 ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 778 ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr); 779 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 780 ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 781 } 782 if (fas->upsmooth) { 783 ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 784 } 785 786 if (fas->downsmooth) { 787 ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 788 } 789 790 if (fas->level != fas->levels - 1) { 791 ierr = SNESSetTolerances(snes, 0.0, 0.0, 0.0, fas->n_cycles, 1000);CHKERRQ(ierr); 792 } 793 794 if (monflg) { 795 fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 796 /* set the monitors for the upsmoother and downsmoother */ 797 ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 798 ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 799 if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 800 if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 801 } else { 802 /* unset the monitors on the coarse levels */ 803 if (fas->level != fas->levels - 1) { 804 ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 805 } 806 } 807 808 /* recursive option setting for the smoothers */ 809 if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 810 PetscFunctionReturn(0); 811 } 812 813 #undef __FUNCT__ 814 #define __FUNCT__ "SNESView_FAS" 815 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 816 { 817 SNES_FAS *fas = (SNES_FAS *) snes->data; 818 PetscBool iascii; 819 PetscErrorCode ierr; 820 PetscInt levels = fas->levels; 821 PetscInt i; 822 823 PetscFunctionBegin; 824 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 825 if (iascii) { 826 ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 827 ierr = PetscViewerASCIIPushTab(viewer); 828 for (i = levels - 1; i >= 0; i--) { 829 ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 830 if (fas->upsmooth) { 831 ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 832 ierr = PetscViewerASCIIPushTab(viewer); 833 ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 834 ierr = PetscViewerASCIIPopTab(viewer); 835 } else { 836 ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 837 } 838 if (fas->downsmooth) { 839 ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 840 ierr = PetscViewerASCIIPushTab(viewer); 841 ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 842 ierr = PetscViewerASCIIPopTab(viewer); 843 } else { 844 ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 845 } 846 if (snes->usegs) { 847 ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n", fas->level);CHKERRQ(ierr); 848 } 849 if (fas->next) fas = (SNES_FAS *)fas->next->data; 850 } 851 ierr = PetscViewerASCIIPopTab(viewer); 852 } else { 853 SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 854 } 855 PetscFunctionReturn(0); 856 } 857 858 /* 859 860 Defines the FAS cycle as: 861 862 fine problem: F(x) = 0 863 coarse problem: F^c(x) = b^c 864 865 b^c = F^c(I^c_fx^f - I^c_fF(x)) 866 867 correction: 868 869 x = x + I(x^c - Rx) 870 871 */ 872 873 #undef __FUNCT__ 874 #define __FUNCT__ "FASCycle_Private" 875 PetscErrorCode FASCycle_Private(SNES snes, Vec X) { 876 877 PetscErrorCode ierr; 878 Vec X_c, Xo_c, F_c, B_c,F,B; 879 SNES_FAS *fas = (SNES_FAS *)snes->data; 880 PetscInt k; 881 PetscReal fnorm; 882 SNESConvergedReason reason; 883 884 PetscFunctionBegin; 885 F = snes->vec_func; 886 B = snes->vec_rhs; 887 /* pre-smooth -- just update using the pre-smoother */ 888 if (fas->downsmooth) { 889 ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 890 /* check convergence reason for the smoother */ 891 ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 892 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 893 snes->reason = SNES_DIVERGED_INNER; 894 PetscFunctionReturn(0); 895 } 896 } else if (snes->usegs && snes->ops->computegs) { 897 if (fas->monitor) { 898 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 899 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 900 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 901 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 902 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 903 } 904 for (k = 0; k < fas->max_up_it; k++) { 905 ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 906 if (fas->monitor) { 907 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 908 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 909 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 910 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 911 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 912 } 913 } 914 } else if (snes->pc) { 915 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 916 ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 917 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 918 snes->reason = SNES_DIVERGED_INNER; 919 PetscFunctionReturn(0); 920 } 921 } 922 if (fas->next) { 923 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 924 925 X_c = fas->next->vec_sol; 926 Xo_c = fas->next->work[0]; 927 F_c = fas->next->vec_func; 928 B_c = fas->next->vec_rhs; 929 930 /* inject the solution */ 931 if (fas->inject) { 932 ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 933 } else { 934 ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 935 ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 936 } 937 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 938 939 /* restrict the defect */ 940 ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 941 942 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 943 fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 944 ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 945 946 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 947 948 /* set initial guess of the coarse problem to the projected fine solution */ 949 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 950 951 /* recurse to the next level */ 952 fas->next->vec_rhs = B_c; 953 /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 954 ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 955 ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 956 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 957 snes->reason = SNES_DIVERGED_INNER; 958 PetscFunctionReturn(0); 959 } 960 /* fas->next->vec_rhs = PETSC_NULL; */ 961 962 /* correct as x <- x + I(x^c - Rx)*/ 963 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 964 ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 965 } 966 /* up-smooth -- just update using the up-smoother */ 967 if (fas->level != 0) { 968 if (fas->upsmooth) { 969 ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 970 ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 971 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 972 snes->reason = SNES_DIVERGED_INNER; 973 PetscFunctionReturn(0); 974 } 975 } else if (snes->usegs && snes->ops->computegs) { 976 if (fas->monitor) { 977 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 978 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 979 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 980 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 981 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 982 } 983 for (k = 0; k < fas->max_down_it; k++) { 984 ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 985 if (fas->monitor) { 986 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 987 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 988 ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 989 ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 990 ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 991 } 992 } 993 } else if (snes->pc) { 994 ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 995 ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 996 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 997 snes->reason = SNES_DIVERGED_INNER; 998 PetscFunctionReturn(0); 999 } 1000 } 1001 } 1002 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1003 1004 PetscFunctionReturn(0); 1005 } 1006 1007 #undef __FUNCT__ 1008 #define __FUNCT__ "SNESSolve_FAS" 1009 1010 PetscErrorCode SNESSolve_FAS(SNES snes) 1011 { 1012 PetscErrorCode ierr; 1013 PetscInt i, maxits; 1014 Vec X, B,F; 1015 PetscReal fnorm; 1016 PetscFunctionBegin; 1017 maxits = snes->max_its; /* maximum number of iterations */ 1018 snes->reason = SNES_CONVERGED_ITERATING; 1019 X = snes->vec_sol; 1020 B = snes->vec_rhs; 1021 F = snes->vec_func; 1022 1023 /*norm setup */ 1024 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1025 snes->iter = 0; 1026 snes->norm = 0.; 1027 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1028 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1029 if (snes->domainerror) { 1030 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1031 PetscFunctionReturn(0); 1032 } 1033 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1034 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1035 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1036 snes->norm = fnorm; 1037 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1038 SNESLogConvHistory(snes,fnorm,0); 1039 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1040 1041 /* set parameter for default relative tolerance convergence test */ 1042 snes->ttol = fnorm*snes->rtol; 1043 /* test convergence */ 1044 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1045 if (snes->reason) PetscFunctionReturn(0); 1046 for (i = 0; i < maxits; i++) { 1047 /* Call general purpose update function */ 1048 1049 if (snes->ops->update) { 1050 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1051 } 1052 ierr = FASCycle_Private(snes, X);CHKERRQ(ierr); 1053 1054 /* check for FAS cycle divergence */ 1055 if (snes->reason != SNES_CONVERGED_ITERATING) { 1056 PetscFunctionReturn(0); 1057 } 1058 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1059 /* Monitor convergence */ 1060 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1061 snes->iter = i+1; 1062 snes->norm = fnorm; 1063 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1064 SNESLogConvHistory(snes,snes->norm,0); 1065 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1066 /* Test for convergence */ 1067 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1068 if (snes->reason) break; 1069 } 1070 if (i == maxits) { 1071 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1072 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1073 } 1074 PetscFunctionReturn(0); 1075 } 1076