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