1 #include "../src/snes/impls/fas/fasimpls.h" /*I "petscsnesfas.h" I*/ 2 3 4 extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *); 5 6 /* -------------- functions called on the fine level -------------- */ 7 8 #undef __FUNCT__ 9 #define __FUNCT__ "SNESFASSetType" 10 /*@ 11 SNESFASSetType - Sets the update and correction type used for FAS. 12 13 Logically Collective 14 15 Input Parameters: 16 snes - FAS context 17 fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE 18 19 .seealso: PCMGSetType() 20 @*/ 21 PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 22 { 23 SNES_FAS *fas = (SNES_FAS*)snes->data; 24 PetscErrorCode ierr; 25 PetscFunctionBegin; 26 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 27 PetscValidLogicalCollectiveEnum(snes,fastype,2); 28 fas->fastype = fastype; 29 if (fas->next) ierr = SNESFASSetType(fas->next, fastype);CHKERRQ(ierr); 30 PetscFunctionReturn(0); 31 } 32 33 34 #undef __FUNCT__ 35 #define __FUNCT__ "SNESFASGetType" 36 /*@ 37 SNESFASGetType - Sets the update and correction type used for FAS. 38 39 Logically Collective 40 41 Input Parameters: 42 .snes - FAS context 43 44 Output Parameters: 45 .fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE 46 47 .seealso: PCMGSetType() 48 @*/ 49 PetscErrorCode SNESFASGetType(SNES snes,SNESFASType *fastype) 50 { 51 SNES_FAS *fas = (SNES_FAS*)snes->data; 52 53 PetscFunctionBegin; 54 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 55 PetscValidPointer(fastype, 2); 56 *fastype = fas->fastype; 57 PetscFunctionReturn(0); 58 } 59 60 #undef __FUNCT__ 61 #define __FUNCT__ "SNESFASSetLevels" 62 /*@C 63 SNESFASSetLevels - Sets the number of levels to use with FAS. 64 Must be called before any other FAS routine. 65 66 Input Parameters: 67 + snes - the snes context 68 . levels - the number of levels 69 - comms - optional communicators for each level; this is to allow solving the coarser 70 problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 71 Fortran. 72 73 Level: intermediate 74 75 Notes: 76 If the number of levels is one then the multigrid uses the -fas_levels prefix 77 for setting the level options rather than the -fas_coarse prefix. 78 79 .keywords: FAS, MG, set, levels, multigrid 80 81 .seealso: SNESFASGetLevels() 82 @*/ 83 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 84 PetscErrorCode ierr; 85 PetscInt i; 86 const char *optionsprefix; 87 char tprefix[128]; 88 SNES_FAS * fas = (SNES_FAS *)snes->data; 89 SNES prevsnes; 90 MPI_Comm comm; 91 PetscFunctionBegin; 92 comm = ((PetscObject)snes)->comm; 93 if (levels == fas->levels) { 94 if (!comms) 95 PetscFunctionReturn(0); 96 } 97 /* user has changed the number of levels; reset */ 98 ierr = SNESReset(snes);CHKERRQ(ierr); 99 /* destroy any coarser levels if necessary */ 100 if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 101 fas->next = PETSC_NULL; 102 fas->previous = PETSC_NULL; 103 prevsnes = snes; 104 /* setup the finest level */ 105 ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 106 for (i = levels - 1; i >= 0; i--) { 107 if (comms) comm = comms[i]; 108 fas->level = i; 109 fas->levels = levels; 110 fas->fine = snes; 111 fas->next = PETSC_NULL; 112 if (i > 0) { 113 ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 114 sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level); 115 ierr = SNESAppendOptionsPrefix(fas->next,tprefix);CHKERRQ(ierr); 116 ierr = SNESAppendOptionsPrefix(fas->next,optionsprefix);CHKERRQ(ierr); 117 ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 118 ierr = SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);CHKERRQ(ierr); 119 ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 120 ((SNES_FAS *)fas->next->data)->previous = prevsnes; 121 prevsnes = fas->next; 122 fas = (SNES_FAS *)prevsnes->data; 123 } 124 } 125 PetscFunctionReturn(0); 126 } 127 128 129 #undef __FUNCT__ 130 #define __FUNCT__ "SNESFASGetLevels" 131 /*@ 132 SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 133 134 Input Parameter: 135 . snes - the nonlinear solver context 136 137 Output parameter: 138 . levels - the number of levels 139 140 Level: advanced 141 142 .keywords: MG, get, levels, multigrid 143 144 .seealso: SNESFASSetLevels(), PCMGGetLevels() 145 @*/ 146 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 147 SNES_FAS * fas = (SNES_FAS *)snes->data; 148 PetscFunctionBegin; 149 *levels = fas->levels; 150 PetscFunctionReturn(0); 151 } 152 153 154 #undef __FUNCT__ 155 #define __FUNCT__ "SNESFASGetCycleSNES" 156 /*@ 157 SNESFASGetCycleSNES - Gets the SNES corresponding to a particular 158 level of the FAS hierarchy. 159 160 Input Parameters: 161 + snes - the multigrid context 162 level - the level to get 163 - lsnes - whether to use the nonlinear smoother or not 164 165 Level: advanced 166 167 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 168 169 .seealso: SNESFASSetLevels(), SNESFASGetLevels() 170 @*/ 171 PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes) { 172 SNES_FAS *fas = (SNES_FAS*)snes->data; 173 PetscInt i; 174 175 PetscFunctionBegin; 176 if (level > fas->levels-1) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels); 177 if (fas->level != fas->levels - 1) 178 SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level); 179 180 *lsnes = snes; 181 for (i = fas->level; i > level; i--) { 182 *lsnes = fas->next; 183 fas = (SNES_FAS *)(*lsnes)->data; 184 } 185 if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 186 PetscFunctionReturn(0); 187 } 188 189 #undef __FUNCT__ 190 #define __FUNCT__ "SNESFASSetNumberSmoothUp" 191 /*@ 192 SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 193 use on all levels. 194 195 Logically Collective on SNES 196 197 Input Parameters: 198 + snes - the multigrid context 199 - n - the number of smoothing steps 200 201 Options Database Key: 202 . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 203 204 Level: advanced 205 206 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 207 208 .seealso: SNESFASSetNumberSmoothDown() 209 @*/ 210 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 211 SNES_FAS * fas = (SNES_FAS *)snes->data; 212 PetscErrorCode ierr = 0; 213 PetscFunctionBegin; 214 fas->max_up_it = n; 215 if (!fas->smoothu) { 216 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu); 217 } 218 ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr); 219 if (fas->next) { 220 ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 221 } 222 PetscFunctionReturn(0); 223 } 224 225 #undef __FUNCT__ 226 #define __FUNCT__ "SNESFASSetNumberSmoothDown" 227 /*@ 228 SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 229 use on all levels. 230 231 Logically Collective on SNES 232 233 Input Parameters: 234 + snes - the multigrid context 235 - n - the number of smoothing steps 236 237 Options Database Key: 238 . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 239 240 Level: advanced 241 242 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 243 244 .seealso: SNESFASSetNumberSmoothUp() 245 @*/ 246 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 247 SNES_FAS * fas = (SNES_FAS *)snes->data; 248 PetscErrorCode ierr = 0; 249 PetscFunctionBegin; 250 if (!fas->smoothu) { 251 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu); 252 } 253 if (!fas->smoothd) { 254 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd); 255 } 256 ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr); 257 fas->max_down_it = n; 258 if (fas->next) { 259 ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 260 } 261 PetscFunctionReturn(0); 262 } 263 264 265 #undef __FUNCT__ 266 #define __FUNCT__ "SNESFASSetCycles" 267 /*@ 268 SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 269 complicated cycling. 270 271 Logically Collective on SNES 272 273 Input Parameters: 274 + snes - the multigrid context 275 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 276 277 Options Database Key: 278 $ -snes_fas_cycles 1 or 2 279 280 Level: advanced 281 282 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 283 284 .seealso: SNESFASSetCyclesOnLevel() 285 @*/ 286 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 287 SNES_FAS * fas = (SNES_FAS *)snes->data; 288 PetscErrorCode ierr; 289 PetscBool isFine; 290 PetscFunctionBegin; 291 ierr = SNESFASCycleIsFine(snes, &isFine); 292 fas->n_cycles = cycles; 293 if (!isFine) 294 ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 295 if (fas->next) { 296 ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 297 } 298 PetscFunctionReturn(0); 299 } 300 301 302 #undef __FUNCT__ 303 #define __FUNCT__ "SNESFASSetMonitor" 304 /*@ 305 SNESFASSetMonitor - Sets the method-specific cycle monitoring 306 307 Logically Collective on SNES 308 309 Input Parameters: 310 + snes - the FAS context 311 - flg - monitor or not 312 313 Level: advanced 314 315 .keywords: FAS, monitor 316 317 .seealso: SNESFASSetCyclesOnLevel() 318 @*/ 319 PetscErrorCode SNESFASSetMonitor(SNES snes, PetscBool flg) { 320 SNES_FAS *fas = (SNES_FAS *)snes->data; 321 PetscErrorCode ierr; 322 PetscBool isFine; 323 PetscInt i, levels = fas->levels; 324 SNES levelsnes; 325 PetscFunctionBegin; 326 ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 327 if (isFine) { 328 for (i = 0; i < levels; i++) { 329 ierr = SNESFASGetCycleSNES(snes, i, &levelsnes); 330 fas = (SNES_FAS *)levelsnes->data; 331 if (flg) { 332 fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)levelsnes)->comm);CHKERRQ(ierr); 333 /* set the monitors for the upsmoother and downsmoother */ 334 ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 335 ierr = SNESMonitorSet(levelsnes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 336 } else { 337 /* unset the monitors on the coarse levels */ 338 if (i != fas->levels - 1) { 339 ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 340 } 341 } 342 } 343 } 344 PetscFunctionReturn(0); 345 } 346 347 #undef __FUNCT__ 348 #define __FUNCT__ "SNESFASCycleCreateSmoother_Private" 349 /* 350 Creates the default smoother type. 351 352 This is SNESNRICHARDSON on each fine level and SNESLS on the coarse level. 353 354 */ 355 PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) { 356 SNES_FAS *fas; 357 const char *optionsprefix; 358 char tprefix[128]; 359 PetscErrorCode ierr; 360 SNES nsmooth; 361 PetscFunctionBegin; 362 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 363 fas = (SNES_FAS*)snes->data; 364 ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 365 /* create the default smoother */ 366 ierr = SNESCreate(((PetscObject)snes)->comm, &nsmooth);CHKERRQ(ierr); 367 if (fas->level == 0) { 368 sprintf(tprefix,"fas_coarse_"); 369 ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 370 ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 371 ierr = SNESSetType(nsmooth, SNESLS);CHKERRQ(ierr); 372 ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr); 373 } else { 374 sprintf(tprefix,"fas_levels_%d_",(int)fas->level); 375 ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 376 ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 377 ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr); 378 ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr); 379 } 380 ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 381 *smooth = nsmooth; 382 PetscFunctionReturn(0); 383 } 384 385 /* ------------- Functions called on a particular level ----------------- */ 386 387 #undef __FUNCT__ 388 #define __FUNCT__ "SNESFASCycleSetCycles" 389 /*@ 390 SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 391 392 Logically Collective on SNES 393 394 Input Parameters: 395 + snes - the multigrid context 396 . level - the level to set the number of cycles on 397 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 398 399 Level: advanced 400 401 .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid 402 403 .seealso: SNESFASSetCycles() 404 @*/ 405 PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) { 406 SNES_FAS * fas = (SNES_FAS *)snes->data; 407 PetscErrorCode ierr; 408 PetscFunctionBegin; 409 fas->n_cycles = cycles; 410 ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 411 PetscFunctionReturn(0); 412 } 413 414 415 #undef __FUNCT__ 416 #define __FUNCT__ "SNESFASCycleGetSmoother" 417 /*@ 418 SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 419 420 Logically Collective on SNES 421 422 Input Parameters: 423 . snes - the multigrid context 424 425 Output Parameters: 426 . smooth - the smoother 427 428 Level: advanced 429 430 .keywords: SNES, FAS, get, smoother, multigrid 431 432 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown() 433 @*/ 434 PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 435 { 436 SNES_FAS *fas; 437 PetscFunctionBegin; 438 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 439 fas = (SNES_FAS*)snes->data; 440 *smooth = fas->smoothd; 441 PetscFunctionReturn(0); 442 } 443 #undef __FUNCT__ 444 #define __FUNCT__ "SNESFASCycleGetSmootherUp" 445 /*@ 446 SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 447 448 Logically Collective on SNES 449 450 Input Parameters: 451 . snes - the multigrid context 452 453 Output Parameters: 454 . smoothu - the smoother 455 456 Notes: 457 Returns the downsmoother if no up smoother is available. This enables transparent 458 default behavior in the process of the solve. 459 460 Level: advanced 461 462 .keywords: SNES, FAS, get, smoother, multigrid 463 464 .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown() 465 @*/ 466 PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 467 { 468 SNES_FAS *fas; 469 PetscFunctionBegin; 470 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 471 fas = (SNES_FAS*)snes->data; 472 if (!fas->smoothu) { 473 *smoothu = fas->smoothd; 474 } else { 475 *smoothu = fas->smoothu; 476 } 477 PetscFunctionReturn(0); 478 } 479 480 #undef __FUNCT__ 481 #define __FUNCT__ "SNESFASCycleGetSmootherDown" 482 /*@ 483 SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 484 485 Logically Collective on SNES 486 487 Input Parameters: 488 . snes - the multigrid context 489 490 Output Parameters: 491 . smoothd - the smoother 492 493 Level: advanced 494 495 .keywords: SNES, FAS, get, smoother, multigrid 496 497 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 498 @*/ 499 PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 500 { 501 SNES_FAS *fas; 502 PetscFunctionBegin; 503 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 504 fas = (SNES_FAS*)snes->data; 505 *smoothd = fas->smoothd; 506 PetscFunctionReturn(0); 507 } 508 509 510 #undef __FUNCT__ 511 #define __FUNCT__ "SNESFASCycleGetCorrection" 512 /*@ 513 SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 514 515 Logically Collective on SNES 516 517 Input Parameters: 518 . snes - the multigrid context 519 520 Output Parameters: 521 . correction - the coarse correction on this level 522 523 Notes: 524 Returns PETSC_NULL on the coarsest level. 525 526 Level: advanced 527 528 .keywords: SNES, FAS, get, smoother, multigrid 529 530 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 531 @*/ 532 PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 533 { 534 SNES_FAS *fas; 535 PetscFunctionBegin; 536 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 537 fas = (SNES_FAS*)snes->data; 538 *correction = fas->next; 539 PetscFunctionReturn(0); 540 } 541 542 #undef __FUNCT__ 543 #define __FUNCT__ "SNESFASCycleGetInterpolation" 544 /*@ 545 SNESFASCycleGetInterpolation - Gets the interpolation on this level 546 547 Logically Collective on SNES 548 549 Input Parameters: 550 . snes - the multigrid context 551 552 Output Parameters: 553 . mat - the interpolation operator on this level 554 555 Level: developer 556 557 .keywords: SNES, FAS, get, smoother, multigrid 558 559 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 560 @*/ 561 PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 562 { 563 SNES_FAS *fas; 564 PetscFunctionBegin; 565 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 566 fas = (SNES_FAS*)snes->data; 567 *mat = fas->interpolate; 568 PetscFunctionReturn(0); 569 } 570 571 572 #undef __FUNCT__ 573 #define __FUNCT__ "SNESFASCycleGetRestriction" 574 /*@ 575 SNESFASCycleGetRestriction - Gets the restriction on this level 576 577 Logically Collective on SNES 578 579 Input Parameters: 580 . snes - the multigrid context 581 582 Output Parameters: 583 . mat - the restriction operator on this level 584 585 Level: developer 586 587 .keywords: SNES, FAS, get, smoother, multigrid 588 589 .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation() 590 @*/ 591 PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 592 { 593 SNES_FAS *fas; 594 PetscFunctionBegin; 595 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 596 fas = (SNES_FAS*)snes->data; 597 *mat = fas->restrct; 598 PetscFunctionReturn(0); 599 } 600 601 602 #undef __FUNCT__ 603 #define __FUNCT__ "SNESFASCycleGetInjection" 604 /*@ 605 SNESFASCycleGetInjection - Gets the injection on this level 606 607 Logically Collective on SNES 608 609 Input Parameters: 610 . snes - the multigrid context 611 612 Output Parameters: 613 . mat - the restriction operator on this level 614 615 Level: developer 616 617 .keywords: SNES, FAS, get, smoother, multigrid 618 619 .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction() 620 @*/ 621 PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 622 { 623 SNES_FAS *fas; 624 PetscFunctionBegin; 625 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 626 fas = (SNES_FAS*)snes->data; 627 *mat = fas->inject; 628 PetscFunctionReturn(0); 629 } 630 631 #undef __FUNCT__ 632 #define __FUNCT__ "SNESFASCycleGetRScale" 633 /*@ 634 SNESFASCycleGetRScale - Gets the injection on this level 635 636 Logically Collective on SNES 637 638 Input Parameters: 639 . snes - the multigrid context 640 641 Output Parameters: 642 . mat - the restriction operator on this level 643 644 Level: developer 645 646 .keywords: SNES, FAS, get, smoother, multigrid 647 648 .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale() 649 @*/ 650 PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 651 { 652 SNES_FAS *fas; 653 PetscFunctionBegin; 654 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 655 fas = (SNES_FAS*)snes->data; 656 *vec = fas->rscale; 657 PetscFunctionReturn(0); 658 } 659 660 #undef __FUNCT__ 661 #define __FUNCT__ "SNESFASCycleIsFine" 662 /*@ 663 SNESFASCycleIsFine - Determines if a given cycle is the fine level. 664 665 Logically Collective on SNES 666 667 Input Parameters: 668 . snes - the FAS context 669 670 Output Parameters: 671 . flg - indicates if this is the fine level or not 672 673 Level: advanced 674 675 .keywords: SNES, FAS 676 677 .seealso: SNESFASSetLevels() 678 @*/ 679 PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 680 { 681 SNES_FAS *fas; 682 PetscFunctionBegin; 683 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 684 fas = (SNES_FAS*)snes->data; 685 if (fas->level == fas->levels - 1) { 686 *flg = PETSC_TRUE; 687 } else { 688 *flg = PETSC_FALSE; 689 } 690 PetscFunctionReturn(0); 691 } 692 693 /* ---------- functions called on the finest level that return level-specific information ---------- */ 694 695 #undef __FUNCT__ 696 #define __FUNCT__ "SNESFASSetInterpolation" 697 /*@ 698 SNESFASSetInterpolation - Sets the function to be used to calculate the 699 interpolation from l-1 to the lth level 700 701 Input Parameters: 702 + snes - the multigrid context 703 . mat - the interpolation operator 704 - level - the level (0 is coarsest) to supply [do not supply 0] 705 706 Level: advanced 707 708 Notes: 709 Usually this is the same matrix used also to set the restriction 710 for the same level. 711 712 One can pass in the interpolation matrix or its transpose; PETSc figures 713 out from the matrix size which one it is. 714 715 .keywords: FAS, multigrid, set, interpolate, level 716 717 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 718 @*/ 719 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 720 SNES_FAS *fas = (SNES_FAS *)snes->data; 721 PetscErrorCode ierr; 722 SNES levelsnes; 723 PetscFunctionBegin; 724 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 725 fas = (SNES_FAS *)levelsnes->data; 726 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 727 ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 728 fas->interpolate = mat; 729 PetscFunctionReturn(0); 730 } 731 732 #undef __FUNCT__ 733 #define __FUNCT__ "SNESFASGetInterpolation" 734 /*@ 735 SNESFASGetInterpolation - Gets the matrix used to calculate the 736 interpolation from l-1 to the lth level 737 738 Input Parameters: 739 + snes - the multigrid context 740 - level - the level (0 is coarsest) to supply [do not supply 0] 741 742 Output Parameters: 743 . mat - the interpolation operator 744 745 Level: advanced 746 747 .keywords: FAS, multigrid, get, interpolate, level 748 749 .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale() 750 @*/ 751 PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) { 752 SNES_FAS *fas = (SNES_FAS *)snes->data; 753 PetscErrorCode ierr; 754 SNES levelsnes; 755 PetscFunctionBegin; 756 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 757 fas = (SNES_FAS *)levelsnes->data; 758 *mat = fas->interpolate; 759 PetscFunctionReturn(0); 760 } 761 762 #undef __FUNCT__ 763 #define __FUNCT__ "SNESFASSetRestriction" 764 /*@ 765 SNESFASSetRestriction - Sets the function to be used to restrict the defect 766 from level l to l-1. 767 768 Input Parameters: 769 + snes - the multigrid context 770 . mat - the restriction matrix 771 - level - the level (0 is coarsest) to supply [Do not supply 0] 772 773 Level: advanced 774 775 Notes: 776 Usually this is the same matrix used also to set the interpolation 777 for the same level. 778 779 One can pass in the interpolation matrix or its transpose; PETSc figures 780 out from the matrix size which one it is. 781 782 If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 783 is used. 784 785 .keywords: FAS, MG, set, multigrid, restriction, level 786 787 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 788 @*/ 789 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 790 SNES_FAS * fas = (SNES_FAS *)snes->data; 791 PetscErrorCode ierr; 792 SNES levelsnes; 793 PetscFunctionBegin; 794 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 795 fas = (SNES_FAS *)levelsnes->data; 796 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 797 ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 798 fas->restrct = mat; 799 PetscFunctionReturn(0); 800 } 801 802 #undef __FUNCT__ 803 #define __FUNCT__ "SNESFASGetRestriction" 804 /*@ 805 SNESFASGetRestriction - Gets the matrix used to calculate the 806 restriction from l to the l-1th level 807 808 Input Parameters: 809 + snes - the multigrid context 810 - level - the level (0 is coarsest) to supply [do not supply 0] 811 812 Output Parameters: 813 . mat - the interpolation operator 814 815 Level: advanced 816 817 .keywords: FAS, multigrid, get, restrict, level 818 819 .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale() 820 @*/ 821 PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) { 822 SNES_FAS *fas = (SNES_FAS *)snes->data; 823 PetscErrorCode ierr; 824 SNES levelsnes; 825 PetscFunctionBegin; 826 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 827 fas = (SNES_FAS *)levelsnes->data; 828 *mat = fas->restrct; 829 PetscFunctionReturn(0); 830 } 831 832 833 #undef __FUNCT__ 834 #define __FUNCT__ "SNESFASSetInjection" 835 /*@ 836 SNESFASSetInjection - Sets the function to be used to inject the solution 837 from level l to l-1. 838 839 Input Parameters: 840 + snes - the multigrid context 841 . mat - the restriction matrix 842 - level - the level (0 is coarsest) to supply [Do not supply 0] 843 844 Level: advanced 845 846 Notes: 847 If you do not set this, the restriction and rscale is used to 848 project the solution instead. 849 850 .keywords: FAS, MG, set, multigrid, restriction, level 851 852 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 853 @*/ 854 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 855 SNES_FAS *fas = (SNES_FAS *)snes->data; 856 PetscErrorCode ierr; 857 SNES levelsnes; 858 PetscFunctionBegin; 859 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 860 fas = (SNES_FAS *)levelsnes->data; 861 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 862 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 863 fas->inject = mat; 864 PetscFunctionReturn(0); 865 } 866 867 868 #undef __FUNCT__ 869 #define __FUNCT__ "SNESFASGetInjection" 870 /*@ 871 SNESFASGetInjection - Gets the matrix used to calculate the 872 injection from l-1 to the lth level 873 874 Input Parameters: 875 + snes - the multigrid context 876 - level - the level (0 is coarsest) to supply [do not supply 0] 877 878 Output Parameters: 879 . mat - the injection operator 880 881 Level: advanced 882 883 .keywords: FAS, multigrid, get, restrict, level 884 885 .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale() 886 @*/ 887 PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) { 888 SNES_FAS *fas = (SNES_FAS *)snes->data; 889 PetscErrorCode ierr; 890 SNES levelsnes; 891 PetscFunctionBegin; 892 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 893 fas = (SNES_FAS *)levelsnes->data; 894 *mat = fas->inject; 895 PetscFunctionReturn(0); 896 } 897 898 #undef __FUNCT__ 899 #define __FUNCT__ "SNESFASSetRScale" 900 /*@ 901 SNESFASSetRScale - Sets the scaling factor of the restriction 902 operator from level l to l-1. 903 904 Input Parameters: 905 + snes - the multigrid context 906 . rscale - the restriction scaling 907 - level - the level (0 is coarsest) to supply [Do not supply 0] 908 909 Level: advanced 910 911 Notes: 912 This is only used in the case that the injection is not set. 913 914 .keywords: FAS, MG, set, multigrid, restriction, level 915 916 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 917 @*/ 918 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 919 SNES_FAS * fas; 920 PetscErrorCode ierr; 921 SNES levelsnes; 922 PetscFunctionBegin; 923 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 924 fas = (SNES_FAS *)levelsnes->data; 925 ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 926 ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 927 fas->rscale = rscale; 928 PetscFunctionReturn(0); 929 } 930 931 #undef __FUNCT__ 932 #define __FUNCT__ "SNESFASGetSmoother" 933 /*@ 934 SNESFASGetSmoother - Gets the default smoother on a level. 935 936 Input Parameters: 937 + snes - the multigrid context 938 - level - the level (0 is coarsest) to supply 939 940 Output Parameters: 941 smooth - the smoother 942 943 Level: advanced 944 945 .keywords: FAS, MG, get, multigrid, smoother, level 946 947 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 948 @*/ 949 PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) { 950 SNES_FAS * fas; 951 PetscErrorCode ierr; 952 SNES levelsnes; 953 PetscFunctionBegin; 954 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 955 fas = (SNES_FAS *)levelsnes->data; 956 if (!fas->smoothd) { 957 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 958 } 959 *smooth = fas->smoothd; 960 PetscFunctionReturn(0); 961 } 962 963 #undef __FUNCT__ 964 #define __FUNCT__ "SNESFASGetSmootherDown" 965 /*@ 966 SNESFASGetSmootherDown - Gets the downsmoother on a level. 967 968 Input Parameters: 969 + snes - the multigrid context 970 - level - the level (0 is coarsest) to supply 971 972 Output Parameters: 973 smooth - the smoother 974 975 Level: advanced 976 977 .keywords: FAS, MG, get, multigrid, smoother, level 978 979 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 980 @*/ 981 PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) { 982 SNES_FAS * fas; 983 PetscErrorCode ierr; 984 SNES levelsnes; 985 PetscFunctionBegin; 986 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 987 fas = (SNES_FAS *)levelsnes->data; 988 /* if the user chooses to differentiate smoothers, create them both at this point */ 989 if (!fas->smoothd) { 990 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 991 } 992 if (!fas->smoothu) { 993 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 994 } 995 *smooth = fas->smoothd; 996 PetscFunctionReturn(0); 997 } 998 999 #undef __FUNCT__ 1000 #define __FUNCT__ "SNESFASGetSmootherUp" 1001 /*@ 1002 SNESFASGetSmootherUp - Gets the upsmoother on a level. 1003 1004 Input Parameters: 1005 + snes - the multigrid context 1006 - level - the level (0 is coarsest) 1007 1008 Output Parameters: 1009 smooth - the smoother 1010 1011 Level: advanced 1012 1013 .keywords: FAS, MG, get, multigrid, smoother, level 1014 1015 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1016 @*/ 1017 PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) { 1018 SNES_FAS * fas; 1019 PetscErrorCode ierr; 1020 SNES levelsnes; 1021 PetscFunctionBegin; 1022 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1023 fas = (SNES_FAS *)levelsnes->data; 1024 /* if the user chooses to differentiate smoothers, create them both at this point */ 1025 if (!fas->smoothd) { 1026 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 1027 } 1028 if (!fas->smoothu) { 1029 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 1030 } 1031 *smooth = fas->smoothu; 1032 PetscFunctionReturn(0); 1033 } 1034 1035 #undef __FUNCT__ 1036 #define __FUNCT__ "SNESFASGetCoarseSolve" 1037 /*@ 1038 SNESFASGetCoarseSolve - Gets the coarsest solver. 1039 1040 Input Parameters: 1041 + snes - the multigrid context 1042 1043 Output Parameters: 1044 solve - the coarse-level solver 1045 1046 Level: advanced 1047 1048 .keywords: FAS, MG, get, multigrid, solver, coarse 1049 1050 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1051 @*/ 1052 PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *smooth) { 1053 SNES_FAS * fas; 1054 PetscErrorCode ierr; 1055 SNES levelsnes; 1056 PetscFunctionBegin; 1057 ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr); 1058 fas = (SNES_FAS *)levelsnes->data; 1059 /* if the user chooses to differentiate smoothers, create them both at this point */ 1060 if (!fas->smoothd) { 1061 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 1062 } 1063 *smooth = fas->smoothd; 1064 PetscFunctionReturn(0); 1065 } 1066