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