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