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