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