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