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