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