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