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