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