1 /* Defines the basic SNES object */ 2 #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3 4 const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0}; 5 6 /*MC 7 8 SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 9 10 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 11 12 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 13 M*/ 14 15 extern PetscErrorCode SNESDestroy_FAS(SNES snes); 16 extern PetscErrorCode SNESSetUp_FAS(SNES snes); 17 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 18 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 19 extern PetscErrorCode SNESSolve_FAS(SNES snes); 20 extern PetscErrorCode SNESReset_FAS(SNES snes); 21 extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 22 23 EXTERN_C_BEGIN 24 #undef __FUNCT__ 25 #define __FUNCT__ "SNESLineSearchSetType_FAS" 26 PetscErrorCode SNESLineSearchSetType_FAS(SNES snes, SNESLineSearchType type) 27 { 28 PetscErrorCode ierr; 29 PetscFunctionBegin; 30 31 switch (type) { 32 case SNES_LS_BASIC: 33 ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 34 break; 35 case SNES_LS_BASIC_NONORMS: 36 ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 37 break; 38 case SNES_LS_QUADRATIC: 39 ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 40 break; 41 default: 42 SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type."); 43 break; 44 } 45 snes->ls_type = type; 46 PetscFunctionReturn(0); 47 } 48 EXTERN_C_END 49 50 EXTERN_C_BEGIN 51 52 #undef __FUNCT__ 53 #define __FUNCT__ "SNESCreate_FAS" 54 PetscErrorCode SNESCreate_FAS(SNES snes) 55 { 56 SNES_FAS * fas; 57 PetscErrorCode ierr; 58 59 PetscFunctionBegin; 60 snes->ops->destroy = SNESDestroy_FAS; 61 snes->ops->setup = SNESSetUp_FAS; 62 snes->ops->setfromoptions = SNESSetFromOptions_FAS; 63 snes->ops->view = SNESView_FAS; 64 snes->ops->solve = SNESSolve_FAS; 65 snes->ops->reset = SNESReset_FAS; 66 67 snes->usesksp = PETSC_FALSE; 68 snes->usespc = PETSC_FALSE; 69 70 ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 71 snes->data = (void*) fas; 72 fas->level = 0; 73 fas->levels = 1; 74 fas->n_cycles = 1; 75 fas->max_up_it = 1; 76 fas->max_down_it = 1; 77 fas->upsmooth = PETSC_NULL; 78 fas->downsmooth = PETSC_NULL; 79 fas->next = PETSC_NULL; 80 fas->previous = PETSC_NULL; 81 fas->interpolate = PETSC_NULL; 82 fas->restrct = PETSC_NULL; 83 fas->inject = PETSC_NULL; 84 fas->monitor = PETSC_NULL; 85 fas->usedmfornumberoflevels = PETSC_FALSE; 86 fas->fastype = SNES_FAS_MULTIPLICATIVE; 87 88 ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_FAS",SNESLineSearchSetType_FAS);CHKERRQ(ierr); 89 ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 90 91 PetscFunctionReturn(0); 92 } 93 EXTERN_C_END 94 95 #undef __FUNCT__ 96 #define __FUNCT__ "SNESFASGetLevels" 97 /*@ 98 SNESFASGetLevels - Gets the number of levels in a FAS. 99 100 Input Parameter: 101 . snes - the preconditioner context 102 103 Output parameter: 104 . levels - the number of levels 105 106 Level: advanced 107 108 .keywords: MG, get, levels, multigrid 109 110 .seealso: SNESFASSetLevels(), PCMGGetLevels() 111 @*/ 112 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 113 SNES_FAS * fas = (SNES_FAS *)snes->data; 114 PetscFunctionBegin; 115 *levels = fas->levels; 116 PetscFunctionReturn(0); 117 } 118 119 #undef __FUNCT__ 120 #define __FUNCT__ "SNESFASSetCycles" 121 /*@ 122 SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 123 complicated cycling. 124 125 Logically Collective on SNES 126 127 Input Parameters: 128 + snes - the multigrid context 129 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 130 131 Options Database Key: 132 $ -snes_fas_cycles 1 or 2 133 134 Level: advanced 135 136 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 137 138 .seealso: SNESFASSetCyclesOnLevel() 139 @*/ 140 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 141 SNES_FAS * fas = (SNES_FAS *)snes->data; 142 PetscErrorCode ierr; 143 PetscFunctionBegin; 144 fas->n_cycles = cycles; 145 if (fas->next) { 146 ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 147 } 148 PetscFunctionReturn(0); 149 } 150 151 #undef __FUNCT__ 152 #define __FUNCT__ "SNESFASSetCyclesOnLevel" 153 /*@ 154 SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 155 156 Logically Collective on SNES 157 158 Input Parameters: 159 + snes - the multigrid context 160 . level - the level to set the number of cycles on 161 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 162 163 Level: advanced 164 165 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 166 167 .seealso: SNESFASSetCycles() 168 @*/ 169 PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 170 SNES_FAS * fas = (SNES_FAS *)snes->data; 171 PetscInt top_level = fas->level,i; 172 173 PetscFunctionBegin; 174 if (level > top_level) 175 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 176 /* get to the correct level */ 177 for (i = fas->level; i > level; i--) { 178 fas = (SNES_FAS *)fas->next->data; 179 } 180 if (fas->level != level) 181 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 182 fas->n_cycles = cycles; 183 PetscFunctionReturn(0); 184 } 185 186 #undef __FUNCT__ 187 #define __FUNCT__ "SNESFASSetGS" 188 /*@C 189 SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 190 Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 191 and nonlinear preconditioners. 192 193 Logically Collective on SNES 194 195 Input Parameters: 196 + snes - the multigrid context 197 . gsfunc - the nonlinear smoother function 198 . ctx - the user context for the nonlinear smoother 199 - use_gs - whether to use the nonlinear smoother or not 200 201 Level: advanced 202 203 .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 204 205 .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 206 @*/ 207 PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 208 PetscErrorCode ierr = 0; 209 SNES_FAS *fas = (SNES_FAS *)snes->data; 210 PetscFunctionBegin; 211 212 if (gsfunc) { 213 ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 214 /* push the provided GS up the tree */ 215 if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 216 } else if (snes->ops->computegs) { 217 /* assume that the user has set the GS solver at this level */ 218 if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 219 } else if (use_gs) { 220 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 221 } 222 PetscFunctionReturn(0); 223 } 224 225 #undef __FUNCT__ 226 #define __FUNCT__ "SNESFASSetGSOnLevel" 227 /*@C 228 SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 229 230 Logically Collective on SNES 231 232 Input Parameters: 233 + snes - the multigrid context 234 . level - the level to set the nonlinear smoother on 235 . gsfunc - the nonlinear smoother function 236 . ctx - the user context for the nonlinear smoother 237 - use_gs - whether to use the nonlinear smoother or not 238 239 Level: advanced 240 241 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 242 243 .seealso: SNESSetGS(), SNESFASSetGS() 244 @*/ 245 PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 246 SNES_FAS *fas = (SNES_FAS *)snes->data; 247 PetscErrorCode ierr; 248 PetscInt top_level = fas->level,i; 249 SNES cur_snes = snes; 250 PetscFunctionBegin; 251 if (level > top_level) 252 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 253 /* get to the correct level */ 254 for (i = fas->level; i > level; i--) { 255 fas = (SNES_FAS *)fas->next->data; 256 cur_snes = fas->next; 257 } 258 if (fas->level != level) 259 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 260 snes->usegs = use_gs; 261 if (gsfunc) { 262 ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 263 } 264 PetscFunctionReturn(0); 265 } 266 267 #undef __FUNCT__ 268 #define __FUNCT__ "SNESFASGetSNES" 269 /*@ 270 SNESFASGetSNES - Gets the SNES corresponding to a particular 271 level of the FAS hierarchy. 272 273 Input Parameters: 274 + snes - the multigrid context 275 level - the level to get 276 - lsnes - whether to use the nonlinear smoother or not 277 278 Level: advanced 279 280 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 281 282 .seealso: SNESFASSetLevels(), SNESFASGetLevels() 283 @*/ 284 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 285 SNES_FAS * fas = (SNES_FAS *)snes->data; 286 PetscInt levels = fas->level; 287 PetscInt i; 288 PetscFunctionBegin; 289 *lsnes = snes; 290 if (fas->level < level) { 291 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 292 } 293 if (level > levels - 1) { 294 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 295 } 296 for (i = fas->level; i > level; i--) { 297 *lsnes = fas->next; 298 fas = (SNES_FAS *)(*lsnes)->data; 299 } 300 if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 301 PetscFunctionReturn(0); 302 } 303 304 #undef __FUNCT__ 305 #define __FUNCT__ "SNESFASSetType" 306 /*@ 307 SNESFASSetType - Sets the update and correction type used for FAS. 308 e 309 310 311 @*/ 312 PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 313 { 314 SNES_FAS *fas = (SNES_FAS*)snes->data; 315 316 PetscFunctionBegin; 317 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 318 PetscValidLogicalCollectiveEnum(snes,fastype,2); 319 fas->fastype = fastype; 320 PetscFunctionReturn(0); 321 } 322 323 324 325 #undef __FUNCT__ 326 #define __FUNCT__ "SNESFASSetLevels" 327 /*@C 328 SNESFASSetLevels - Sets the number of levels to use with FAS. 329 Must be called before any other FAS routine. 330 331 Input Parameters: 332 + snes - the snes context 333 . levels - the number of levels 334 - comms - optional communicators for each level; this is to allow solving the coarser 335 problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 336 Fortran. 337 338 Level: intermediate 339 340 Notes: 341 If the number of levels is one then the multigrid uses the -fas_levels prefix 342 for setting the level options rather than the -fas_coarse prefix. 343 344 .keywords: FAS, MG, set, levels, multigrid 345 346 .seealso: SNESFASGetLevels() 347 @*/ 348 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 349 PetscErrorCode ierr; 350 PetscInt i; 351 SNES_FAS * fas = (SNES_FAS *)snes->data; 352 SNES prevsnes; 353 MPI_Comm comm; 354 PetscFunctionBegin; 355 comm = ((PetscObject)snes)->comm; 356 if (levels == fas->levels) { 357 if (!comms) 358 PetscFunctionReturn(0); 359 } 360 /* user has changed the number of levels; reset */ 361 ierr = SNESReset(snes);CHKERRQ(ierr); 362 /* destroy any coarser levels if necessary */ 363 if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 364 fas->next = PETSC_NULL; 365 fas->previous = PETSC_NULL; 366 prevsnes = snes; 367 /* setup the finest level */ 368 for (i = levels - 1; i >= 0; i--) { 369 if (comms) comm = comms[i]; 370 fas->level = i; 371 fas->levels = levels; 372 fas->next = PETSC_NULL; 373 if (i > 0) { 374 ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 375 ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 376 ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 377 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 378 ((SNES_FAS *)fas->next->data)->previous = prevsnes; 379 prevsnes = fas->next; 380 fas = (SNES_FAS *)prevsnes->data; 381 } 382 } 383 PetscFunctionReturn(0); 384 } 385 386 #undef __FUNCT__ 387 #define __FUNCT__ "SNESFASSetNumberSmoothUp" 388 /*@ 389 SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 390 use on all levels. 391 392 Logically Collective on SNES 393 394 Input Parameters: 395 + snes - the multigrid context 396 - n - the number of smoothing steps 397 398 Options Database Key: 399 . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 400 401 Level: advanced 402 403 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 404 405 .seealso: SNESFASSetNumberSmoothDown() 406 @*/ 407 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 408 SNES_FAS * fas = (SNES_FAS *)snes->data; 409 PetscErrorCode ierr = 0; 410 PetscFunctionBegin; 411 fas->max_up_it = n; 412 if (fas->next) { 413 ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 414 } 415 PetscFunctionReturn(0); 416 } 417 418 #undef __FUNCT__ 419 #define __FUNCT__ "SNESFASSetNumberSmoothDown" 420 /*@ 421 SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 422 use on all levels. 423 424 Logically Collective on SNES 425 426 Input Parameters: 427 + snes - the multigrid context 428 - n - the number of smoothing steps 429 430 Options Database Key: 431 . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 432 433 Level: advanced 434 435 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 436 437 .seealso: SNESFASSetNumberSmoothUp() 438 @*/ 439 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 440 SNES_FAS * fas = (SNES_FAS *)snes->data; 441 PetscErrorCode ierr = 0; 442 PetscFunctionBegin; 443 fas->max_down_it = n; 444 if (fas->next) { 445 ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 446 } 447 PetscFunctionReturn(0); 448 } 449 450 #undef __FUNCT__ 451 #define __FUNCT__ "SNESFASSetInterpolation" 452 /*@ 453 SNESFASSetInterpolation - Sets the function to be used to calculate the 454 interpolation from l-1 to the lth level 455 456 Input Parameters: 457 + snes - the multigrid context 458 . mat - the interpolation operator 459 - level - the level (0 is coarsest) to supply [do not supply 0] 460 461 Level: advanced 462 463 Notes: 464 Usually this is the same matrix used also to set the restriction 465 for the same level. 466 467 One can pass in the interpolation matrix or its transpose; PETSc figures 468 out from the matrix size which one it is. 469 470 .keywords: FAS, multigrid, set, interpolate, level 471 472 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale() 473 @*/ 474 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 475 SNES_FAS * fas = (SNES_FAS *)snes->data; 476 PetscInt top_level = fas->level,i; 477 478 PetscFunctionBegin; 479 if (level > top_level) 480 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 481 /* get to the correct level */ 482 for (i = fas->level; i > level; i--) { 483 fas = (SNES_FAS *)fas->next->data; 484 } 485 if (fas->level != level) 486 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 487 fas->interpolate = mat; 488 PetscFunctionReturn(0); 489 } 490 491 #undef __FUNCT__ 492 #define __FUNCT__ "SNESFASSetRestriction" 493 /*@ 494 SNESFASSetRestriction - Sets the function to be used to restrict the defect 495 from level l to l-1. 496 497 Input Parameters: 498 + snes - the multigrid context 499 . mat - the restriction matrix 500 - level - the level (0 is coarsest) to supply [Do not supply 0] 501 502 Level: advanced 503 504 Notes: 505 Usually this is the same matrix used also to set the interpolation 506 for the same level. 507 508 One can pass in the interpolation matrix or its transpose; PETSc figures 509 out from the matrix size which one it is. 510 511 If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 512 is used. 513 514 .keywords: FAS, MG, set, multigrid, restriction, level 515 516 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 517 @*/ 518 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 519 SNES_FAS * fas = (SNES_FAS *)snes->data; 520 PetscInt top_level = fas->level,i; 521 522 PetscFunctionBegin; 523 if (level > top_level) 524 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 525 /* get to the correct level */ 526 for (i = fas->level; i > level; i--) { 527 fas = (SNES_FAS *)fas->next->data; 528 } 529 if (fas->level != level) 530 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 531 fas->restrct = mat; 532 PetscFunctionReturn(0); 533 } 534 535 536 #undef __FUNCT__ 537 #define __FUNCT__ "SNESFASSetRScale" 538 /*@ 539 SNESFASSetRScale - Sets the scaling factor of the restriction 540 operator from level l to l-1. 541 542 Input Parameters: 543 + snes - the multigrid context 544 . rscale - the restriction scaling 545 - level - the level (0 is coarsest) to supply [Do not supply 0] 546 547 Level: advanced 548 549 Notes: 550 This is only used in the case that the injection is not set. 551 552 .keywords: FAS, MG, set, multigrid, restriction, level 553 554 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 555 @*/ 556 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 557 SNES_FAS * fas = (SNES_FAS *)snes->data; 558 PetscInt top_level = fas->level,i; 559 560 PetscFunctionBegin; 561 if (level > top_level) 562 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 563 /* get to the correct level */ 564 for (i = fas->level; i > level; i--) { 565 fas = (SNES_FAS *)fas->next->data; 566 } 567 if (fas->level != level) 568 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 569 fas->rscale = rscale; 570 PetscFunctionReturn(0); 571 } 572 573 574 #undef __FUNCT__ 575 #define __FUNCT__ "SNESFASSetInjection" 576 /*@ 577 SNESFASSetInjection - Sets the function to be used to inject the solution 578 from level l to l-1. 579 580 Input Parameters: 581 + snes - the multigrid context 582 . mat - the restriction matrix 583 - level - the level (0 is coarsest) to supply [Do not supply 0] 584 585 Level: advanced 586 587 Notes: 588 If you do not set this, the restriction and rscale is used to 589 project the solution instead. 590 591 .keywords: FAS, MG, set, multigrid, restriction, level 592 593 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 594 @*/ 595 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 596 SNES_FAS * fas = (SNES_FAS *)snes->data; 597 PetscInt top_level = fas->level,i; 598 599 PetscFunctionBegin; 600 if (level > top_level) 601 SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 602 /* get to the correct level */ 603 for (i = fas->level; i > level; i--) { 604 fas = (SNES_FAS *)fas->next->data; 605 } 606 if (fas->level != level) 607 SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 608 fas->inject = mat; 609 PetscFunctionReturn(0); 610 } 611 612 #undef __FUNCT__ 613 #define __FUNCT__ "SNESReset_FAS" 614 PetscErrorCode SNESReset_FAS(SNES snes) 615 { 616 PetscErrorCode ierr = 0; 617 SNES_FAS * fas = (SNES_FAS *)snes->data; 618 619 PetscFunctionBegin; 620 if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 621 if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 622 if (fas->inject) { 623 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 624 } 625 if (fas->interpolate == fas->restrct) { 626 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 627 fas->restrct = PETSC_NULL; 628 } else { 629 if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 630 if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 631 } 632 if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 633 634 if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 635 if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 636 if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 637 PetscFunctionReturn(0); 638 } 639 640 #undef __FUNCT__ 641 #define __FUNCT__ "SNESDestroy_FAS" 642 PetscErrorCode SNESDestroy_FAS(SNES snes) 643 { 644 SNES_FAS * fas = (SNES_FAS *)snes->data; 645 PetscErrorCode ierr = 0; 646 647 PetscFunctionBegin; 648 /* recursively resets and then destroys */ 649 if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 650 if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 651 if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 652 ierr = PetscFree(fas);CHKERRQ(ierr); 653 654 PetscFunctionReturn(0); 655 } 656 657 #undef __FUNCT__ 658 #define __FUNCT__ "SNESSetUp_FAS" 659 PetscErrorCode SNESSetUp_FAS(SNES snes) 660 { 661 SNES_FAS *fas = (SNES_FAS *) snes->data; 662 PetscErrorCode ierr; 663 VecScatter injscatter; 664 PetscInt dm_levels; 665 Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 666 667 PetscFunctionBegin; 668 669 if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 670 ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 671 dm_levels++; 672 if (dm_levels > fas->levels) { 673 674 /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESSetUp_FAS*/ 675 vec_sol = snes->vec_sol; 676 vec_func = snes->vec_func; 677 vec_sol_update = snes->vec_sol_update; 678 vec_rhs = snes->vec_rhs; 679 snes->vec_sol = PETSC_NULL; 680 snes->vec_func = PETSC_NULL; 681 snes->vec_sol_update = PETSC_NULL; 682 snes->vec_rhs = PETSC_NULL; 683 684 /* reset the number of levels */ 685 ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 686 ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 687 688 snes->vec_sol = vec_sol; 689 snes->vec_func = vec_func; 690 snes->vec_rhs = vec_rhs; 691 snes->vec_sol_update = vec_sol_update; 692 } 693 } 694 695 if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 696 697 if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 698 ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 699 } else { 700 ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 701 } 702 703 /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 704 if (fas->upsmooth) { 705 if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 706 ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 707 } 708 if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 709 ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 710 } 711 if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 712 ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 713 } 714 ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 715 } 716 if (fas->downsmooth) { 717 if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 718 ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 719 } 720 if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 721 ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 722 } 723 if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 724 ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 725 } 726 ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 727 } 728 /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 729 if (fas->next) { 730 if (fas->galerkin) { 731 ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 732 } else { 733 if (snes->ops->computefunction && !fas->next->ops->computefunction) { 734 ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 735 } 736 if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 737 ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 738 } 739 if (snes->ops->computegs && !fas->next->ops->computegs) { 740 ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 741 } 742 } 743 } 744 if (snes->dm) { 745 /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 746 if (fas->next) { 747 /* for now -- assume the DM and the evaluation functions have been set externally */ 748 if (!fas->next->dm) { 749 ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 750 ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 751 } 752 /* set the interpolation and restriction from the DM */ 753 if (!fas->interpolate) { 754 ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 755 fas->restrct = fas->interpolate; 756 } 757 /* set the injection from the DM */ 758 if (!fas->inject) { 759 ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 760 ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 761 ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 762 } 763 } 764 /* set the DMs of the pre and post-smoothers here */ 765 if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 766 if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 767 } 768 769 /* setup FAS work vectors */ 770 if (fas->galerkin) { 771 ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 772 ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 773 } 774 775 if (fas->next) { 776 /* gotta set up the solution vector for this to work */ 777 if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 778 if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);} 779 ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 780 } 781 /* got to set them all up at once */ 782 PetscFunctionReturn(0); 783 } 784 785 #undef __FUNCT__ 786 #define __FUNCT__ "SNESSetFromOptions_FAS" 787 PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 788 { 789 SNES_FAS *fas = (SNES_FAS *) snes->data; 790 PetscInt levels = 1; 791 PetscBool flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg; 792 PetscErrorCode ierr; 793 const char *def_smooth = SNESNRICHARDSON; 794 char pre_type[256]; 795 char post_type[256]; 796 char monfilename[PETSC_MAX_PATH_LEN]; 797 SNESFASType fastype; 798 const char *optionsprefix; 799 const char *prefix; 800 801 PetscFunctionBegin; 802 ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 803 804 /* number of levels -- only process on the finest level */ 805 if (fas->levels - 1 == fas->level) { 806 ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 807 if (!flg && snes->dm) { 808 ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 809 levels++; 810 fas->usedmfornumberoflevels = PETSC_TRUE; 811 } 812 ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 813 } 814 815 /* type of pre and/or post smoothers -- set both at once */ 816 ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 817 ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 818 fastype = fas->fastype; 819 ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 820 if (flg) { 821 ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 822 } 823 824 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 825 826 /* smoother setup options */ 827 ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr); 828 ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr); 829 ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr); 830 if (fas->level == 0) { 831 ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr); 832 } 833 /* options for the number of preconditioning cycles and cycle type */ 834 835 ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 836 837 ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 838 839 ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 840 841 842 ierr = PetscOptionsTail();CHKERRQ(ierr); 843 /* setup from the determined types if there is no pointwise procedure or smoother defined */ 844 845 if ((!fas->downsmooth) && smoothcoarseflg) { 846 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 847 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 848 ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 849 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 850 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 851 } 852 853 if ((!fas->downsmooth) && smoothdownflg) { 854 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 855 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 856 ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 857 if (fas->level || (fas->levels == 1)) { 858 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr); 859 } else { 860 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 861 } 862 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 863 } 864 865 if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) { 866 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 867 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 868 ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 869 ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr); 870 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 871 } 872 873 if ((!fas->downsmooth) && smoothflg) { 874 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 875 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 876 ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 877 if (fas->level || (fas->levels == 1)) { 878 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr); 879 } else { 880 ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 881 } 882 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 883 } 884 885 if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) { 886 ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 887 ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 888 ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 889 ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr); 890 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 891 } 892 893 if (fas->upsmooth) { 894 ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr); 895 } 896 897 if (fas->downsmooth) { 898 ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr); 899 } 900 901 if (fas->level != fas->levels - 1) { 902 ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr); 903 } 904 905 if (monflg) { 906 fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 907 /* set the monitors for the upsmoother and downsmoother */ 908 ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 909 ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 910 if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 911 if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 912 } else { 913 /* unset the monitors on the coarse levels */ 914 if (fas->level != fas->levels - 1) { 915 ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 916 } 917 } 918 919 /* recursive option setting for the smoothers */ 920 if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 921 PetscFunctionReturn(0); 922 } 923 924 #undef __FUNCT__ 925 #define __FUNCT__ "SNESView_FAS" 926 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 927 { 928 SNES_FAS *fas = (SNES_FAS *) snes->data; 929 PetscBool iascii; 930 PetscErrorCode ierr; 931 932 PetscFunctionBegin; 933 ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 934 if (iascii) { 935 ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 936 ierr = PetscViewerASCIIPushTab(viewer); 937 ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 938 if (fas->upsmooth) { 939 ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 940 ierr = PetscViewerASCIIPushTab(viewer); 941 ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 942 ierr = PetscViewerASCIIPopTab(viewer); 943 } else { 944 ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 945 } 946 if (fas->downsmooth) { 947 ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 948 ierr = PetscViewerASCIIPushTab(viewer); 949 ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 950 ierr = PetscViewerASCIIPopTab(viewer); 951 } else { 952 ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 953 } 954 ierr = PetscViewerASCIIPopTab(viewer); 955 } else { 956 SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 957 } 958 PetscFunctionReturn(0); 959 } 960 961 #undef __FUNCT__ 962 #define __FUNCT__ "FASDownSmooth" 963 /* 964 Defines the action of the downsmoother 965 */ 966 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 967 PetscErrorCode ierr = 0; 968 PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 969 SNESConvergedReason reason; 970 SNES_FAS *fas = (SNES_FAS *)snes->data; 971 Vec G, W, Y; 972 PetscBool lssuccess; 973 PetscInt k; 974 PetscFunctionBegin; 975 G = snes->work[1]; 976 W = snes->work[2]; 977 Y = snes->work[3]; 978 if (fas->downsmooth) { 979 ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 980 /* check convergence reason for the smoother */ 981 ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 982 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 983 snes->reason = SNES_DIVERGED_INNER; 984 PetscFunctionReturn(0); 985 } 986 } else { 987 /* basic richardson smoother */ 988 for (k = 0; k < fas->max_up_it; k++) { 989 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 990 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 991 ierr = VecCopy(F, Y);CHKERRQ(ierr); 992 ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 993 if (!lssuccess) { 994 if (++snes->numFailures >= snes->maxFailures) { 995 snes->reason = SNES_DIVERGED_LINE_SEARCH; 996 PetscFunctionReturn(0); 997 } 998 } 999 ierr = VecCopy(W, X);CHKERRQ(ierr); 1000 ierr = VecCopy(G, F);CHKERRQ(ierr); 1001 fnorm = gnorm; 1002 } 1003 } 1004 PetscFunctionReturn(0); 1005 } 1006 1007 1008 #undef __FUNCT__ 1009 #define __FUNCT__ "FASUpSmooth" 1010 /* 1011 Defines the action of the upsmoother 1012 */ 1013 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 1014 PetscErrorCode ierr = 0; 1015 PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 1016 SNESConvergedReason reason; 1017 SNES_FAS *fas = (SNES_FAS *)snes->data; 1018 Vec G, W, Y; 1019 PetscBool lssuccess; 1020 PetscInt k; 1021 PetscFunctionBegin; 1022 G = snes->work[1]; 1023 W = snes->work[2]; 1024 Y = snes->work[3]; 1025 if (fas->upsmooth) { 1026 ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 1027 /* check convergence reason for the smoother */ 1028 ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 1029 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1030 snes->reason = SNES_DIVERGED_INNER; 1031 PetscFunctionReturn(0); 1032 } 1033 } else { 1034 /* basic richardson smoother */ 1035 for (k = 0; k < fas->max_up_it; k++) { 1036 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1037 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1038 ierr = VecCopy(F, Y);CHKERRQ(ierr); 1039 ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 1040 if (!lssuccess) { 1041 if (++snes->numFailures >= snes->maxFailures) { 1042 snes->reason = SNES_DIVERGED_LINE_SEARCH; 1043 PetscFunctionReturn(0); 1044 } 1045 } 1046 ierr = VecCopy(W, X);CHKERRQ(ierr); 1047 ierr = VecCopy(G, F);CHKERRQ(ierr); 1048 fnorm = gnorm; 1049 } 1050 } 1051 PetscFunctionReturn(0); 1052 } 1053 1054 #undef __FUNCT__ 1055 #define __FUNCT__ "FASCoarseCorrection" 1056 /* 1057 1058 Performs the FAS coarse correction as: 1059 1060 fine problem: F(x) = 0 1061 coarse problem: F^c(x) = b^c 1062 1063 b^c = F^c(I^c_fx^f - I^c_fF(x)) 1064 1065 with correction: 1066 1067 1068 1069 */ 1070 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 1071 PetscErrorCode ierr; 1072 Vec X_c, Xo_c, F_c, B_c; 1073 SNES_FAS *fas = (SNES_FAS *)snes->data; 1074 SNESConvergedReason reason; 1075 PetscFunctionBegin; 1076 if (fas->next) { 1077 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1078 1079 X_c = fas->next->vec_sol; 1080 Xo_c = fas->next->work[0]; 1081 F_c = fas->next->vec_func; 1082 B_c = fas->next->vec_rhs; 1083 1084 /* inject the solution */ 1085 if (fas->inject) { 1086 ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 1087 } else { 1088 ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 1089 ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 1090 } 1091 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1092 1093 /* restrict the defect */ 1094 ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1095 1096 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1097 fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1098 ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1099 1100 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1101 1102 /* set initial guess of the coarse problem to the projected fine solution */ 1103 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1104 1105 /* recurse to the next level */ 1106 fas->next->vec_rhs = B_c; 1107 /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1108 ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1109 ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1110 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1111 snes->reason = SNES_DIVERGED_INNER; 1112 PetscFunctionReturn(0); 1113 } 1114 /* fas->next->vec_rhs = PETSC_NULL; */ 1115 1116 /* correct as x <- x + I(x^c - Rx)*/ 1117 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1118 ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1119 } 1120 PetscFunctionReturn(0); 1121 } 1122 1123 #undef __FUNCT__ 1124 #define __FUNCT__ "FASCycle_Additive" 1125 /* 1126 1127 The additive cycle looks like: 1128 1129 xhat = x 1130 xhat = dS(x, b) 1131 x = coarsecorrection(xhat, b_d) 1132 x = x + nu*(xhat - x); 1133 (optional) x = uS(x, b) 1134 1135 With the coarse RHS (defect correction) as below. 1136 1137 */ 1138 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 1139 Vec F, B, Xhat; 1140 Vec X_c, Xo_c, F_c, B_c, G, W; 1141 PetscErrorCode ierr; 1142 SNES_FAS * fas = (SNES_FAS *)snes->data; 1143 SNESConvergedReason reason; 1144 PetscReal xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.; 1145 PetscBool lssucceed; 1146 PetscFunctionBegin; 1147 1148 F = snes->vec_func; 1149 B = snes->vec_rhs; 1150 Xhat = snes->work[3]; 1151 G = snes->work[1]; 1152 W = snes->work[2]; 1153 ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 1154 /* recurse first */ 1155 if (fas->next) { 1156 ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 1157 1158 X_c = fas->next->vec_sol; 1159 Xo_c = fas->next->work[0]; 1160 F_c = fas->next->vec_func; 1161 B_c = fas->next->vec_rhs; 1162 1163 /* inject the solution */ 1164 if (fas->inject) { 1165 ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr); 1166 } else { 1167 ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr); 1168 ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 1169 } 1170 ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1171 1172 /* restrict the defect */ 1173 ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1174 1175 /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1176 fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1177 ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1178 1179 ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1180 1181 /* set initial guess of the coarse problem to the projected fine solution */ 1182 ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1183 1184 /* recurse */ 1185 fas->next->vec_rhs = B_c; 1186 ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1187 1188 /* smooth on this level */ 1189 ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 1190 1191 ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1192 if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1193 snes->reason = SNES_DIVERGED_INNER; 1194 PetscFunctionReturn(0); 1195 } 1196 1197 /* correct as x <- x + I(x^c - Rx)*/ 1198 ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1199 ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 1200 1201 /* additive correction of the coarse direction*/ 1202 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1203 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1204 ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr); 1205 ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1206 ierr = VecCopy(W, X);CHKERRQ(ierr); 1207 ierr = VecCopy(G, F);CHKERRQ(ierr); 1208 fnorm = gnorm; 1209 } else { 1210 ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 1211 } 1212 PetscFunctionReturn(0); 1213 } 1214 1215 #undef __FUNCT__ 1216 #define __FUNCT__ "FASCycle_Multiplicative" 1217 /* 1218 1219 Defines the FAS cycle as: 1220 1221 fine problem: F(x) = 0 1222 coarse problem: F^c(x) = b^c 1223 1224 b^c = F^c(I^c_fx^f - I^c_fF(x)) 1225 1226 correction: 1227 1228 x = x + I(x^c - Rx) 1229 1230 */ 1231 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 1232 1233 PetscErrorCode ierr; 1234 Vec F,B; 1235 SNES_FAS *fas = (SNES_FAS *)snes->data; 1236 1237 PetscFunctionBegin; 1238 F = snes->vec_func; 1239 B = snes->vec_rhs; 1240 /* pre-smooth -- just update using the pre-smoother */ 1241 ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 1242 1243 ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 1244 1245 if (fas->level != 0) { 1246 ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1247 } 1248 ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1249 1250 PetscFunctionReturn(0); 1251 } 1252 1253 #undef __FUNCT__ 1254 #define __FUNCT__ "SNESSolve_FAS" 1255 1256 PetscErrorCode SNESSolve_FAS(SNES snes) 1257 { 1258 PetscErrorCode ierr; 1259 PetscInt i, maxits; 1260 Vec X, F; 1261 PetscReal fnorm; 1262 SNES_FAS *fas = (SNES_FAS *)snes->data; 1263 PetscFunctionBegin; 1264 maxits = snes->max_its; /* maximum number of iterations */ 1265 snes->reason = SNES_CONVERGED_ITERATING; 1266 X = snes->vec_sol; 1267 F = snes->vec_func; 1268 1269 /*norm setup */ 1270 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1271 snes->iter = 0; 1272 snes->norm = 0.; 1273 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1274 ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1275 if (snes->domainerror) { 1276 snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1277 PetscFunctionReturn(0); 1278 } 1279 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1280 if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1281 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1282 snes->norm = fnorm; 1283 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1284 SNESLogConvHistory(snes,fnorm,0); 1285 ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1286 1287 /* set parameter for default relative tolerance convergence test */ 1288 snes->ttol = fnorm*snes->rtol; 1289 /* test convergence */ 1290 ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1291 if (snes->reason) PetscFunctionReturn(0); 1292 for (i = 0; i < maxits; i++) { 1293 /* Call general purpose update function */ 1294 1295 if (snes->ops->update) { 1296 ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1297 } 1298 if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 1299 ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 1300 } else { 1301 ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 1302 } 1303 1304 /* check for FAS cycle divergence */ 1305 if (snes->reason != SNES_CONVERGED_ITERATING) { 1306 PetscFunctionReturn(0); 1307 } 1308 ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1309 /* Monitor convergence */ 1310 ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1311 snes->iter = i+1; 1312 snes->norm = fnorm; 1313 ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1314 SNESLogConvHistory(snes,snes->norm,0); 1315 ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1316 /* Test for convergence */ 1317 ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1318 if (snes->reason) break; 1319 } 1320 if (i == maxits) { 1321 ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1322 if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1323 } 1324 PetscFunctionReturn(0); 1325 } 1326