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