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 #undef __FUNCT__ 302 #define __FUNCT__ "SNESFASCycleCreateSmoother_Private" 303 /* 304 Creates the default smoother type. 305 306 This is SNESNRICHARDSON on each fine level and SNESLS on the coarse level. 307 308 */ 309 PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) { 310 SNES_FAS *fas; 311 const char *optionsprefix; 312 char tprefix[128]; 313 PetscErrorCode ierr; 314 SNES nsmooth; 315 PetscFunctionBegin; 316 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 317 fas = (SNES_FAS*)snes->data; 318 ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 319 /* create the default smoother */ 320 ierr = SNESCreate(((PetscObject)snes)->comm, &nsmooth);CHKERRQ(ierr); 321 if (fas->level == 0) { 322 sprintf(tprefix,"fas_coarse_"); 323 ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 324 ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 325 ierr = SNESSetType(nsmooth, SNESLS);CHKERRQ(ierr); 326 ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, fas->max_up_it, 1000);CHKERRQ(ierr); 327 ierr = SNESSetFromOptions(nsmooth);CHKERRQ(ierr); 328 } else { 329 sprintf(tprefix,"fas_levels_%d_",(int)fas->level); 330 ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 331 ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 332 ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr); 333 ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, fas->max_up_it, 1000);CHKERRQ(ierr); 334 ierr = SNESSetFromOptions(nsmooth);CHKERRQ(ierr); 335 } 336 ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 337 *smooth = nsmooth; 338 PetscFunctionReturn(0); 339 } 340 341 /* ------------- Functions called on a particular level ----------------- */ 342 343 #undef __FUNCT__ 344 #define __FUNCT__ "SNESFASCycleSetCycles" 345 /*@ 346 SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 347 348 Logically Collective on SNES 349 350 Input Parameters: 351 + snes - the multigrid context 352 . level - the level to set the number of cycles on 353 - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 354 355 Level: advanced 356 357 .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid 358 359 .seealso: SNESFASSetCycles() 360 @*/ 361 PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) { 362 SNES_FAS * fas = (SNES_FAS *)snes->data; 363 PetscErrorCode ierr; 364 PetscFunctionBegin; 365 fas->n_cycles = cycles; 366 ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 367 PetscFunctionReturn(0); 368 } 369 370 371 #undef __FUNCT__ 372 #define __FUNCT__ "SNESFASCycleGetSmoother" 373 /*@ 374 SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 375 376 Logically Collective on SNES 377 378 Input Parameters: 379 . snes - the multigrid context 380 381 Output Parameters: 382 . smooth - the smoother 383 384 Level: advanced 385 386 .keywords: SNES, FAS, get, smoother, multigrid 387 388 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown() 389 @*/ 390 PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 391 { 392 SNES_FAS *fas; 393 PetscFunctionBegin; 394 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 395 fas = (SNES_FAS*)snes->data; 396 *smooth = fas->smoothd; 397 PetscFunctionReturn(0); 398 } 399 #undef __FUNCT__ 400 #define __FUNCT__ "SNESFASCycleGetSmootherUp" 401 /*@ 402 SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 403 404 Logically Collective on SNES 405 406 Input Parameters: 407 . snes - the multigrid context 408 409 Output Parameters: 410 . smoothu - the smoother 411 412 Notes: 413 Returns the downsmoother if no up smoother is available. This enables transparent 414 default behavior in the process of the solve. 415 416 Level: advanced 417 418 .keywords: SNES, FAS, get, smoother, multigrid 419 420 .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown() 421 @*/ 422 PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 423 { 424 SNES_FAS *fas; 425 PetscFunctionBegin; 426 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 427 fas = (SNES_FAS*)snes->data; 428 if (!fas->smoothu) { 429 *smoothu = fas->smoothd; 430 } else { 431 *smoothu = fas->smoothu; 432 } 433 PetscFunctionReturn(0); 434 } 435 436 #undef __FUNCT__ 437 #define __FUNCT__ "SNESFASCycleGetSmootherDown" 438 /*@ 439 SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 440 441 Logically Collective on SNES 442 443 Input Parameters: 444 . snes - the multigrid context 445 446 Output Parameters: 447 . smoothd - the smoother 448 449 Level: advanced 450 451 .keywords: SNES, FAS, get, smoother, multigrid 452 453 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 454 @*/ 455 PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 456 { 457 SNES_FAS *fas; 458 PetscFunctionBegin; 459 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 460 fas = (SNES_FAS*)snes->data; 461 *smoothd = fas->smoothd; 462 PetscFunctionReturn(0); 463 } 464 465 466 #undef __FUNCT__ 467 #define __FUNCT__ "SNESFASCycleGetCorrection" 468 /*@ 469 SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 470 471 Logically Collective on SNES 472 473 Input Parameters: 474 . snes - the multigrid context 475 476 Output Parameters: 477 . correction - the coarse correction on this level 478 479 Notes: 480 Returns PETSC_NULL on the coarsest level. 481 482 Level: advanced 483 484 .keywords: SNES, FAS, get, smoother, multigrid 485 486 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 487 @*/ 488 PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 489 { 490 SNES_FAS *fas; 491 PetscFunctionBegin; 492 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 493 fas = (SNES_FAS*)snes->data; 494 *correction = fas->next; 495 PetscFunctionReturn(0); 496 } 497 498 #undef __FUNCT__ 499 #define __FUNCT__ "SNESFASCycleGetInterpolation" 500 /*@ 501 SNESFASCycleGetInterpolation - Gets the interpolation on this level 502 503 Logically Collective on SNES 504 505 Input Parameters: 506 . snes - the multigrid context 507 508 Output Parameters: 509 . mat - the interpolation operator on this level 510 511 Level: developer 512 513 .keywords: SNES, FAS, get, smoother, multigrid 514 515 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 516 @*/ 517 PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 518 { 519 SNES_FAS *fas; 520 PetscFunctionBegin; 521 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 522 fas = (SNES_FAS*)snes->data; 523 *mat = fas->interpolate; 524 PetscFunctionReturn(0); 525 } 526 527 528 #undef __FUNCT__ 529 #define __FUNCT__ "SNESFASCycleGetRestriction" 530 /*@ 531 SNESFASCycleGetRestriction - Gets the restriction on this level 532 533 Logically Collective on SNES 534 535 Input Parameters: 536 . snes - the multigrid context 537 538 Output Parameters: 539 . mat - the restriction operator on this level 540 541 Level: developer 542 543 .keywords: SNES, FAS, get, smoother, multigrid 544 545 .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation() 546 @*/ 547 PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 548 { 549 SNES_FAS *fas; 550 PetscFunctionBegin; 551 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 552 fas = (SNES_FAS*)snes->data; 553 *mat = fas->restrct; 554 PetscFunctionReturn(0); 555 } 556 557 558 #undef __FUNCT__ 559 #define __FUNCT__ "SNESFASCycleGetInjection" 560 /*@ 561 SNESFASCycleGetInjection - Gets the injection on this level 562 563 Logically Collective on SNES 564 565 Input Parameters: 566 . snes - the multigrid context 567 568 Output Parameters: 569 . mat - the restriction operator on this level 570 571 Level: developer 572 573 .keywords: SNES, FAS, get, smoother, multigrid 574 575 .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction() 576 @*/ 577 PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 578 { 579 SNES_FAS *fas; 580 PetscFunctionBegin; 581 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 582 fas = (SNES_FAS*)snes->data; 583 *mat = fas->inject; 584 PetscFunctionReturn(0); 585 } 586 587 #undef __FUNCT__ 588 #define __FUNCT__ "SNESFASCycleGetRScale" 589 /*@ 590 SNESFASCycleGetRScale - Gets the injection on this level 591 592 Logically Collective on SNES 593 594 Input Parameters: 595 . snes - the multigrid context 596 597 Output Parameters: 598 . mat - the restriction operator on this level 599 600 Level: developer 601 602 .keywords: SNES, FAS, get, smoother, multigrid 603 604 .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale() 605 @*/ 606 PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 607 { 608 SNES_FAS *fas; 609 PetscFunctionBegin; 610 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 611 fas = (SNES_FAS*)snes->data; 612 *vec = fas->rscale; 613 PetscFunctionReturn(0); 614 } 615 616 #undef __FUNCT__ 617 #define __FUNCT__ "SNESFASCycleIsFine" 618 /*@ 619 SNESFASCycleIsFine - Determines if a given cycle is the fine level. 620 621 Logically Collective on SNES 622 623 Input Parameters: 624 . snes - the FAS context 625 626 Output Parameters: 627 . flg - indicates if this is the fine level or not 628 629 Level: advanced 630 631 .keywords: SNES, FAS 632 633 .seealso: SNESFASSetLevels() 634 @*/ 635 PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 636 { 637 SNES_FAS *fas; 638 PetscFunctionBegin; 639 PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 640 fas = (SNES_FAS*)snes->data; 641 if (fas->level == fas->levels - 1) { 642 *flg = PETSC_TRUE; 643 } else { 644 *flg = PETSC_FALSE; 645 } 646 PetscFunctionReturn(0); 647 } 648 649 /* ---------- functions called on the finest level that return level-specific information ---------- */ 650 651 #undef __FUNCT__ 652 #define __FUNCT__ "SNESFASSetInterpolation" 653 /*@ 654 SNESFASSetInterpolation - Sets the function to be used to calculate the 655 interpolation from l-1 to the lth level 656 657 Input Parameters: 658 + snes - the multigrid context 659 . mat - the interpolation operator 660 - level - the level (0 is coarsest) to supply [do not supply 0] 661 662 Level: advanced 663 664 Notes: 665 Usually this is the same matrix used also to set the restriction 666 for the same level. 667 668 One can pass in the interpolation matrix or its transpose; PETSc figures 669 out from the matrix size which one it is. 670 671 .keywords: FAS, multigrid, set, interpolate, level 672 673 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 674 @*/ 675 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 676 SNES_FAS *fas = (SNES_FAS *)snes->data; 677 PetscErrorCode ierr; 678 SNES levelsnes; 679 PetscFunctionBegin; 680 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 681 fas = (SNES_FAS *)levelsnes->data; 682 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 683 ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 684 fas->interpolate = mat; 685 PetscFunctionReturn(0); 686 } 687 688 #undef __FUNCT__ 689 #define __FUNCT__ "SNESFASGetInterpolation" 690 /*@ 691 SNESFASGetInterpolation - Gets the matrix used to calculate the 692 interpolation from l-1 to the lth level 693 694 Input Parameters: 695 + snes - the multigrid context 696 - level - the level (0 is coarsest) to supply [do not supply 0] 697 698 Output Parameters: 699 . mat - the interpolation operator 700 701 Level: advanced 702 703 .keywords: FAS, multigrid, get, interpolate, level 704 705 .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale() 706 @*/ 707 PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) { 708 SNES_FAS *fas = (SNES_FAS *)snes->data; 709 PetscErrorCode ierr; 710 SNES levelsnes; 711 PetscFunctionBegin; 712 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 713 fas = (SNES_FAS *)levelsnes->data; 714 *mat = fas->interpolate; 715 PetscFunctionReturn(0); 716 } 717 718 #undef __FUNCT__ 719 #define __FUNCT__ "SNESFASSetRestriction" 720 /*@ 721 SNESFASSetRestriction - Sets the function to be used to restrict the defect 722 from level l to l-1. 723 724 Input Parameters: 725 + snes - the multigrid context 726 . mat - the restriction matrix 727 - level - the level (0 is coarsest) to supply [Do not supply 0] 728 729 Level: advanced 730 731 Notes: 732 Usually this is the same matrix used also to set the interpolation 733 for the same level. 734 735 One can pass in the interpolation matrix or its transpose; PETSc figures 736 out from the matrix size which one it is. 737 738 If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 739 is used. 740 741 .keywords: FAS, MG, set, multigrid, restriction, level 742 743 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 744 @*/ 745 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 746 SNES_FAS * fas = (SNES_FAS *)snes->data; 747 PetscErrorCode ierr; 748 SNES levelsnes; 749 PetscFunctionBegin; 750 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 751 fas = (SNES_FAS *)levelsnes->data; 752 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 753 ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 754 fas->restrct = mat; 755 PetscFunctionReturn(0); 756 } 757 758 #undef __FUNCT__ 759 #define __FUNCT__ "SNESFASGetRestriction" 760 /*@ 761 SNESFASGetRestriction - Gets the matrix used to calculate the 762 restriction from l to the l-1th level 763 764 Input Parameters: 765 + snes - the multigrid context 766 - level - the level (0 is coarsest) to supply [do not supply 0] 767 768 Output Parameters: 769 . mat - the interpolation operator 770 771 Level: advanced 772 773 .keywords: FAS, multigrid, get, restrict, level 774 775 .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale() 776 @*/ 777 PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) { 778 SNES_FAS *fas = (SNES_FAS *)snes->data; 779 PetscErrorCode ierr; 780 SNES levelsnes; 781 PetscFunctionBegin; 782 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 783 fas = (SNES_FAS *)levelsnes->data; 784 *mat = fas->restrct; 785 PetscFunctionReturn(0); 786 } 787 788 789 #undef __FUNCT__ 790 #define __FUNCT__ "SNESFASSetInjection" 791 /*@ 792 SNESFASSetInjection - Sets the function to be used to inject the solution 793 from level l to l-1. 794 795 Input Parameters: 796 + snes - the multigrid context 797 . mat - the restriction matrix 798 - level - the level (0 is coarsest) to supply [Do not supply 0] 799 800 Level: advanced 801 802 Notes: 803 If you do not set this, the restriction and rscale is used to 804 project the solution instead. 805 806 .keywords: FAS, MG, set, multigrid, restriction, level 807 808 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 809 @*/ 810 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 811 SNES_FAS *fas = (SNES_FAS *)snes->data; 812 PetscErrorCode ierr; 813 SNES levelsnes; 814 PetscFunctionBegin; 815 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 816 fas = (SNES_FAS *)levelsnes->data; 817 ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 818 ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 819 fas->inject = mat; 820 PetscFunctionReturn(0); 821 } 822 823 824 #undef __FUNCT__ 825 #define __FUNCT__ "SNESFASGetInjection" 826 /*@ 827 SNESFASGetInjection - Gets the matrix used to calculate the 828 injection from l-1 to the lth level 829 830 Input Parameters: 831 + snes - the multigrid context 832 - level - the level (0 is coarsest) to supply [do not supply 0] 833 834 Output Parameters: 835 . mat - the injection operator 836 837 Level: advanced 838 839 .keywords: FAS, multigrid, get, restrict, level 840 841 .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale() 842 @*/ 843 PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) { 844 SNES_FAS *fas = (SNES_FAS *)snes->data; 845 PetscErrorCode ierr; 846 SNES levelsnes; 847 PetscFunctionBegin; 848 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 849 fas = (SNES_FAS *)levelsnes->data; 850 *mat = fas->inject; 851 PetscFunctionReturn(0); 852 } 853 854 #undef __FUNCT__ 855 #define __FUNCT__ "SNESFASSetRScale" 856 /*@ 857 SNESFASSetRScale - Sets the scaling factor of the restriction 858 operator from level l to l-1. 859 860 Input Parameters: 861 + snes - the multigrid context 862 . rscale - the restriction scaling 863 - level - the level (0 is coarsest) to supply [Do not supply 0] 864 865 Level: advanced 866 867 Notes: 868 This is only used in the case that the injection is not set. 869 870 .keywords: FAS, MG, set, multigrid, restriction, level 871 872 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 873 @*/ 874 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 875 SNES_FAS * fas; 876 PetscErrorCode ierr; 877 SNES levelsnes; 878 PetscFunctionBegin; 879 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 880 fas = (SNES_FAS *)levelsnes->data; 881 ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 882 ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 883 fas->rscale = rscale; 884 PetscFunctionReturn(0); 885 } 886 887 #undef __FUNCT__ 888 #define __FUNCT__ "SNESFASGetSmoother" 889 /*@ 890 SNESFASGetSmoother - Gets the default smoother on a level. 891 892 Input Parameters: 893 + snes - the multigrid context 894 - level - the level (0 is coarsest) to supply 895 896 Output Parameters: 897 smooth - the smoother 898 899 Level: advanced 900 901 .keywords: FAS, MG, get, multigrid, smoother, level 902 903 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 904 @*/ 905 PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) { 906 SNES_FAS * fas; 907 PetscErrorCode ierr; 908 SNES levelsnes; 909 PetscFunctionBegin; 910 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 911 fas = (SNES_FAS *)levelsnes->data; 912 if (!fas->smoothd) { 913 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 914 } 915 *smooth = fas->smoothd; 916 PetscFunctionReturn(0); 917 } 918 919 #undef __FUNCT__ 920 #define __FUNCT__ "SNESFASGetSmootherDown" 921 /*@ 922 SNESFASGetSmootherDown - Gets the downsmoother on a level. 923 924 Input Parameters: 925 + snes - the multigrid context 926 - level - the level (0 is coarsest) to supply 927 928 Output Parameters: 929 smooth - the smoother 930 931 Level: advanced 932 933 .keywords: FAS, MG, get, multigrid, smoother, level 934 935 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 936 @*/ 937 PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) { 938 SNES_FAS * fas; 939 PetscErrorCode ierr; 940 SNES levelsnes; 941 PetscFunctionBegin; 942 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 943 fas = (SNES_FAS *)levelsnes->data; 944 /* if the user chooses to differentiate smoothers, create them both at this point */ 945 if (!fas->smoothd) { 946 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 947 } 948 if (!fas->smoothu) { 949 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 950 } 951 *smooth = fas->smoothd; 952 PetscFunctionReturn(0); 953 } 954 955 #undef __FUNCT__ 956 #define __FUNCT__ "SNESFASGetSmootherUp" 957 /*@ 958 SNESFASGetSmootherUp - Gets the upsmoother on a level. 959 960 Input Parameters: 961 + snes - the multigrid context 962 - level - the level (0 is coarsest) 963 964 Output Parameters: 965 smooth - the smoother 966 967 Level: advanced 968 969 .keywords: FAS, MG, get, multigrid, smoother, level 970 971 .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 972 @*/ 973 PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) { 974 SNES_FAS * fas; 975 PetscErrorCode ierr; 976 SNES levelsnes; 977 PetscFunctionBegin; 978 ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 979 fas = (SNES_FAS *)levelsnes->data; 980 /* if the user chooses to differentiate smoothers, create them both at this point */ 981 if (!fas->smoothd) { 982 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 983 } 984 if (!fas->smoothu) { 985 ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 986 } 987 *smooth = fas->smoothu; 988 PetscFunctionReturn(0); 989 } 990