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