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