1a7b5fb5fSBarry Smith #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnes.h" I*/ 2ab8d36c9SPeter Brune 3ab8d36c9SPeter Brune /*@ 4ab8d36c9SPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 5ab8d36c9SPeter Brune 6ab8d36c9SPeter Brune Logically Collective 7ab8d36c9SPeter Brune 8ab8d36c9SPeter Brune Input Parameters: 9583a1250SSatish Balay + snes - FAS context 10*f6dfbefdSBarry Smith - fastype - `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, or `SNES_FAS_KASKADE` 11583a1250SSatish Balay 12583a1250SSatish Balay Level: intermediate 13ab8d36c9SPeter Brune 14*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `PCMGSetType()`, `SNESFASGetType()` 15ab8d36c9SPeter Brune @*/ 169371c9d4SSatish Balay PetscErrorCode SNESFASSetType(SNES snes, SNESFASType fastype) { 17f833ba53SLisandro Dalcin SNES_FAS *fas; 1822d28d08SBarry Smith 19ab8d36c9SPeter Brune PetscFunctionBegin; 20f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 21ab8d36c9SPeter Brune PetscValidLogicalCollectiveEnum(snes, fastype, 2); 22f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 23ab8d36c9SPeter Brune fas->fastype = fastype; 241baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetType(fas->next, fastype)); 25ab8d36c9SPeter Brune PetscFunctionReturn(0); 26ab8d36c9SPeter Brune } 27ab8d36c9SPeter Brune 28ab8d36c9SPeter Brune /*@ 29*f6dfbefdSBarry Smith SNESFASGetType - Gets the update and correction type used for FAS. 30ab8d36c9SPeter Brune 31ab8d36c9SPeter Brune Logically Collective 32ab8d36c9SPeter Brune 33*f6dfbefdSBarry Smith Input Parameter: 34*f6dfbefdSBarry Smith . snes - `SNESFAS` context 35ab8d36c9SPeter Brune 36*f6dfbefdSBarry Smith Output Parameter: 37*f6dfbefdSBarry Smith . fastype - `SNES_FAS_ADDITIVE` or `SNES_FAS_MULTIPLICATIVE` 38ab8d36c9SPeter Brune 39583a1250SSatish Balay Level: intermediate 40583a1250SSatish Balay 41*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `PCMGSetType()`, `SNESFASSetType()` 42ab8d36c9SPeter Brune @*/ 439371c9d4SSatish Balay PetscErrorCode SNESFASGetType(SNES snes, SNESFASType *fastype) { 44f833ba53SLisandro Dalcin SNES_FAS *fas; 45ab8d36c9SPeter Brune 46ab8d36c9SPeter Brune PetscFunctionBegin; 47f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 48ab8d36c9SPeter Brune PetscValidPointer(fastype, 2); 49f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 50ab8d36c9SPeter Brune *fastype = fas->fastype; 51ab8d36c9SPeter Brune PetscFunctionReturn(0); 52ab8d36c9SPeter Brune } 53ab8d36c9SPeter Brune 54ab8d36c9SPeter Brune /*@C 55*f6dfbefdSBarry Smith SNESFASSetLevels - Sets the number of levels to use with `SNESFAS`. 56ab8d36c9SPeter Brune Must be called before any other FAS routine. 57ab8d36c9SPeter Brune 58ab8d36c9SPeter Brune Input Parameters: 59ab8d36c9SPeter Brune + snes - the snes context 60ab8d36c9SPeter Brune . levels - the number of levels 61ab8d36c9SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 622bf68e3eSBarry Smith problems on smaller sets of processors. 63ab8d36c9SPeter Brune 64ab8d36c9SPeter Brune Level: intermediate 65ab8d36c9SPeter Brune 66*f6dfbefdSBarry Smith Note: 67ab8d36c9SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 68ab8d36c9SPeter Brune for setting the level options rather than the -fas_coarse prefix. 69ab8d36c9SPeter Brune 70*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetLevels()` 71ab8d36c9SPeter Brune @*/ 729371c9d4SSatish Balay PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms) { 73ab8d36c9SPeter Brune PetscInt i; 74ab8d36c9SPeter Brune const char *optionsprefix; 75ab8d36c9SPeter Brune char tprefix[128]; 76f833ba53SLisandro Dalcin SNES_FAS *fas; 77ab8d36c9SPeter Brune SNES prevsnes; 78ab8d36c9SPeter Brune MPI_Comm comm; 7922d28d08SBarry Smith 80ab8d36c9SPeter Brune PetscFunctionBegin; 81f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 82f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 839566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 84ab8d36c9SPeter Brune if (levels == fas->levels) { 8522d28d08SBarry Smith if (!comms) PetscFunctionReturn(0); 86ab8d36c9SPeter Brune } 87ab8d36c9SPeter Brune /* user has changed the number of levels; reset */ 88dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, reset); 89ab8d36c9SPeter Brune /* destroy any coarser levels if necessary */ 909566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&fas->next)); 910298fd71SBarry Smith fas->next = NULL; 920298fd71SBarry Smith fas->previous = NULL; 93ab8d36c9SPeter Brune prevsnes = snes; 94ab8d36c9SPeter Brune /* setup the finest level */ 959566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 969566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)snes, PetscMGLevelId, levels - 1)); 97ab8d36c9SPeter Brune for (i = levels - 1; i >= 0; i--) { 98ab8d36c9SPeter Brune if (comms) comm = comms[i]; 99ab8d36c9SPeter Brune fas->level = i; 100ab8d36c9SPeter Brune fas->levels = levels; 101ab8d36c9SPeter Brune fas->fine = snes; 1020298fd71SBarry Smith fas->next = NULL; 103ab8d36c9SPeter Brune if (i > 0) { 1049566063dSJacob Faibussowitsch PetscCall(SNESCreate(comm, &fas->next)); 1059566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 1069566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%d_cycle_", (int)fas->level)); 1079566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->next, optionsprefix)); 1089566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->next, tprefix)); 1099566063dSJacob Faibussowitsch PetscCall(SNESSetType(fas->next, SNESFAS)); 1109566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs)); 1119566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i)); 1129566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)fas->next, PetscMGLevelId, i - 1)); 1131aa26658SKarl Rupp 114ab8d36c9SPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 1151aa26658SKarl Rupp 116ab8d36c9SPeter Brune prevsnes = fas->next; 117ab8d36c9SPeter Brune fas = (SNES_FAS *)prevsnes->data; 118ab8d36c9SPeter Brune } 119ab8d36c9SPeter Brune } 120ab8d36c9SPeter Brune PetscFunctionReturn(0); 121ab8d36c9SPeter Brune } 122ab8d36c9SPeter Brune 123ab8d36c9SPeter Brune /*@ 124*f6dfbefdSBarry Smith SNESFASGetLevels - Gets the number of levels in a `SNESFAS`, including fine and coarse grids 125ab8d36c9SPeter Brune 126ab8d36c9SPeter Brune Input Parameter: 127*f6dfbefdSBarry Smith . snes - the `SNES` nonlinear solver context 128ab8d36c9SPeter Brune 129ab8d36c9SPeter Brune Output parameter: 130ab8d36c9SPeter Brune . levels - the number of levels 131ab8d36c9SPeter Brune 132ab8d36c9SPeter Brune Level: advanced 133ab8d36c9SPeter Brune 134*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()`, `PCMGGetLevels()` 135ab8d36c9SPeter Brune @*/ 1369371c9d4SSatish Balay PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels) { 137f833ba53SLisandro Dalcin SNES_FAS *fas; 1385fd66863SKarl Rupp 139ab8d36c9SPeter Brune PetscFunctionBegin; 140f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 141064a246eSJacob Faibussowitsch PetscValidIntPointer(levels, 2); 142f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 143ab8d36c9SPeter Brune *levels = fas->levels; 144ab8d36c9SPeter Brune PetscFunctionReturn(0); 145ab8d36c9SPeter Brune } 146ab8d36c9SPeter Brune 147ab8d36c9SPeter Brune /*@ 148*f6dfbefdSBarry Smith SNESFASGetCycleSNES - Gets the `SNES` corresponding to a particular 149*f6dfbefdSBarry Smith level of the `SNESFAS` hierarchy. 150ab8d36c9SPeter Brune 151ab8d36c9SPeter Brune Input Parameters: 152*f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context 153*f6dfbefdSBarry Smith - level - the level to get 154*f6dfbefdSBarry Smith 155*f6dfbefdSBarry Smith Output Parameter: 156*f6dfbefdSBarry Smith . lsnes - the `SNES` for the requested level 157ab8d36c9SPeter Brune 158ab8d36c9SPeter Brune Level: advanced 159ab8d36c9SPeter Brune 160*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()`, `SNESFASGetLevels()` 161ab8d36c9SPeter Brune @*/ 1629371c9d4SSatish Balay PetscErrorCode SNESFASGetCycleSNES(SNES snes, PetscInt level, SNES *lsnes) { 163f833ba53SLisandro Dalcin SNES_FAS *fas; 164ab8d36c9SPeter Brune PetscInt i; 165ab8d36c9SPeter Brune 166ab8d36c9SPeter Brune PetscFunctionBegin; 167f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 168f833ba53SLisandro Dalcin PetscValidPointer(lsnes, 3); 169f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 17063a3b9bcSJacob Faibussowitsch 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); 1710b121fc5SBarry Smith PetscCheck(fas->level == fas->levels - 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetCycleSNES may only be called on the finest-level SNES"); 172ab8d36c9SPeter Brune 173ab8d36c9SPeter Brune *lsnes = snes; 174ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) { 175ab8d36c9SPeter Brune *lsnes = fas->next; 176ab8d36c9SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 177ab8d36c9SPeter Brune } 17808401ef6SPierre Jolivet PetscCheck(fas->level == level, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "SNESFAS level hierarchy corrupt"); 179ab8d36c9SPeter Brune PetscFunctionReturn(0); 180ab8d36c9SPeter Brune } 181ab8d36c9SPeter Brune 182ab8d36c9SPeter Brune /*@ 183ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 184ab8d36c9SPeter Brune use on all levels. 185ab8d36c9SPeter Brune 186*f6dfbefdSBarry Smith Logically Collective on snes 187ab8d36c9SPeter Brune 188ab8d36c9SPeter Brune Input Parameters: 189*f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context 190ab8d36c9SPeter Brune - n - the number of smoothing steps 191ab8d36c9SPeter Brune 192ab8d36c9SPeter Brune Options Database Key: 193ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 194ab8d36c9SPeter Brune 195ab8d36c9SPeter Brune Level: advanced 196ab8d36c9SPeter Brune 197*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothDown()` 198ab8d36c9SPeter Brune @*/ 1999371c9d4SSatish Balay PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 200f833ba53SLisandro Dalcin SNES_FAS *fas; 20122d28d08SBarry Smith 202ab8d36c9SPeter Brune PetscFunctionBegin; 203f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 204f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 205ab8d36c9SPeter Brune fas->max_up_it = n; 20648a46eb9SPierre Jolivet if (!fas->smoothu && fas->level != 0) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu)); 2071baa6e33SBarry Smith if (fas->smoothu) PetscCall(SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs)); 2081baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothUp(fas->next, n)); 209ab8d36c9SPeter Brune PetscFunctionReturn(0); 210ab8d36c9SPeter Brune } 211ab8d36c9SPeter Brune 212ab8d36c9SPeter Brune /*@ 213ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 214ab8d36c9SPeter Brune use on all levels. 215ab8d36c9SPeter Brune 216*f6dfbefdSBarry Smith Logically Collective on snes 217ab8d36c9SPeter Brune 218ab8d36c9SPeter Brune Input Parameters: 219*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 220ab8d36c9SPeter Brune - n - the number of smoothing steps 221ab8d36c9SPeter Brune 222ab8d36c9SPeter Brune Options Database Key: 223ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 224ab8d36c9SPeter Brune 225ab8d36c9SPeter Brune Level: advanced 226ab8d36c9SPeter Brune 227*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 228ab8d36c9SPeter Brune @*/ 2299371c9d4SSatish Balay PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 230f833ba53SLisandro Dalcin SNES_FAS *fas; 23122d28d08SBarry Smith 232ab8d36c9SPeter Brune PetscFunctionBegin; 233f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 234f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 23548a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd)); 2369566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs)); 2371aa26658SKarl Rupp 238ab8d36c9SPeter Brune fas->max_down_it = n; 2391baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothDown(fas->next, n)); 240ab8d36c9SPeter Brune PetscFunctionReturn(0); 241ab8d36c9SPeter Brune } 242ab8d36c9SPeter Brune 24387f44e3fSPeter Brune /*@ 244*f6dfbefdSBarry Smith SNESFASSetContinuation - Sets the `SNESFAS` cycle to default to exact Newton solves on the upsweep 24587f44e3fSPeter Brune 246*f6dfbefdSBarry Smith Logically Collective on snes 24787f44e3fSPeter Brune 24887f44e3fSPeter Brune Input Parameters: 249*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 25087f44e3fSPeter Brune - n - the number of smoothing steps 25187f44e3fSPeter Brune 25287f44e3fSPeter Brune Options Database Key: 25387f44e3fSPeter Brune . -snes_fas_continuation - sets continuation to true 25487f44e3fSPeter Brune 25587f44e3fSPeter Brune Level: advanced 25687f44e3fSPeter Brune 257*f6dfbefdSBarry Smith Note: 25895452b02SPatrick Sanan This sets the prefix on the upsweep smoothers to -fas_continuation 25987f44e3fSPeter Brune 260*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 26187f44e3fSPeter Brune @*/ 2629371c9d4SSatish Balay PetscErrorCode SNESFASSetContinuation(SNES snes, PetscBool continuation) { 26387f44e3fSPeter Brune const char *optionsprefix; 26487f44e3fSPeter Brune char tprefix[128]; 265f833ba53SLisandro Dalcin SNES_FAS *fas; 26687f44e3fSPeter Brune 26787f44e3fSPeter Brune PetscFunctionBegin; 268f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 269f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 2709566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 27148a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu)); 2729566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_levels_continuation_", sizeof(tprefix))); 2739566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(fas->smoothu, optionsprefix)); 2749566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->smoothu, tprefix)); 2759566063dSJacob Faibussowitsch PetscCall(SNESSetType(fas->smoothu, SNESNEWTONLS)); 2769566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothu, fas->fine->abstol, fas->fine->rtol, fas->fine->stol, 50, 100)); 27787f44e3fSPeter Brune fas->continuation = continuation; 2781baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetContinuation(fas->next, continuation)); 27987f44e3fSPeter Brune PetscFunctionReturn(0); 28087f44e3fSPeter Brune } 28187f44e3fSPeter Brune 282ab8d36c9SPeter Brune /*@ 283*f6dfbefdSBarry Smith SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use `SNESFASSetCyclesOnLevel()` for more 284ab8d36c9SPeter Brune complicated cycling. 285ab8d36c9SPeter Brune 286*f6dfbefdSBarry Smith Logically Collective on snes 287ab8d36c9SPeter Brune 288ab8d36c9SPeter Brune Input Parameters: 289*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 290ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 291ab8d36c9SPeter Brune 292ab8d36c9SPeter Brune Options Database Key: 29367b8a455SSatish Balay . -snes_fas_cycles <1,2> - 1 for V-cycle, 2 for W-cycle 294ab8d36c9SPeter Brune 295ab8d36c9SPeter Brune Level: advanced 296ab8d36c9SPeter Brune 297*f6dfbefdSBarry Smith .seealso: `SNES`, `SNESFAS`, `SNESFASSetCyclesOnLevel()` 298ab8d36c9SPeter Brune @*/ 2999371c9d4SSatish Balay PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 300f833ba53SLisandro Dalcin SNES_FAS *fas; 301ab8d36c9SPeter Brune PetscBool isFine; 30222d28d08SBarry Smith 303ab8d36c9SPeter Brune PetscFunctionBegin; 304f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3059566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 306f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 307ab8d36c9SPeter Brune fas->n_cycles = cycles; 30848a46eb9SPierre Jolivet if (!isFine) PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs)); 3091baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetCycles(fas->next, cycles)); 310ab8d36c9SPeter Brune PetscFunctionReturn(0); 311ab8d36c9SPeter Brune } 312ab8d36c9SPeter Brune 313c8c899caSPeter Brune /*@ 314c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring 315c8c899caSPeter Brune 316*f6dfbefdSBarry Smith Logically Collective on snes 317c8c899caSPeter Brune 318c8c899caSPeter Brune Input Parameters: 319*f6dfbefdSBarry Smith + snes - the `SNESFAS` context 320d142ab34SLawrence Mitchell . vf - viewer and format structure (may be NULL if flg is FALSE) 321c8c899caSPeter Brune - flg - monitor or not 322c8c899caSPeter Brune 323c8c899caSPeter Brune Level: advanced 324c8c899caSPeter Brune 325*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESSetMonitor()`, `SNESFASSetCyclesOnLevel()` 326c8c899caSPeter Brune @*/ 3279371c9d4SSatish Balay PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg) { 328f833ba53SLisandro Dalcin SNES_FAS *fas; 329c8c899caSPeter Brune PetscBool isFine; 330f833ba53SLisandro Dalcin PetscInt i, levels; 331c8c899caSPeter Brune SNES levelsnes; 33222d28d08SBarry Smith 333c8c899caSPeter Brune PetscFunctionBegin; 334f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3359566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 336f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 337f833ba53SLisandro Dalcin levels = fas->levels; 338c8c899caSPeter Brune if (isFine) { 339c8c899caSPeter Brune for (i = 0; i < levels; i++) { 3409566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 341c8c899caSPeter Brune fas = (SNES_FAS *)levelsnes->data; 342c8c899caSPeter Brune if (flg) { 343c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */ 3449566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes)); 345d142ab34SLawrence Mitchell /* Only register destroy on finest level */ 3469566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(levelsnes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))SNESMonitorDefault, vf, (!i ? (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy : NULL))); 3471aa26658SKarl Rupp } else if (i != fas->levels - 1) { 348c8c899caSPeter Brune /* unset the monitors on the coarse levels */ 3499566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes)); 350c8c899caSPeter Brune } 351c8c899caSPeter Brune } 352c8c899caSPeter Brune } 353c8c899caSPeter Brune PetscFunctionReturn(0); 354c8c899caSPeter Brune } 355c8c899caSPeter Brune 3560dd27c6cSPeter Brune /*@ 357*f6dfbefdSBarry Smith SNESFASSetLog - Sets or unsets time logging for various `SNESFAS` stages on all levels 3580dd27c6cSPeter Brune 359*f6dfbefdSBarry Smith Logically Collective on snes 3600dd27c6cSPeter Brune 3610dd27c6cSPeter Brune Input Parameters: 362*f6dfbefdSBarry Smith + snes - the `SNESFAS` context 3630dd27c6cSPeter Brune - flg - monitor or not 3640dd27c6cSPeter Brune 3650dd27c6cSPeter Brune Level: advanced 3660dd27c6cSPeter Brune 367*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetMonitor()` 3680dd27c6cSPeter Brune @*/ 3699371c9d4SSatish Balay PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) { 370f833ba53SLisandro Dalcin SNES_FAS *fas; 3710dd27c6cSPeter Brune PetscBool isFine; 372f833ba53SLisandro Dalcin PetscInt i, levels; 3730dd27c6cSPeter Brune SNES levelsnes; 3740dd27c6cSPeter Brune char eventname[128]; 3750dd27c6cSPeter Brune 3760dd27c6cSPeter Brune PetscFunctionBegin; 377f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3789566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 379f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 380f833ba53SLisandro Dalcin levels = fas->levels; 3810dd27c6cSPeter Brune if (isFine) { 3820dd27c6cSPeter Brune for (i = 0; i < levels; i++) { 3839566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 3840dd27c6cSPeter Brune fas = (SNES_FAS *)levelsnes->data; 3850dd27c6cSPeter Brune if (flg) { 3869566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSetup %d", (int)i)); 3879566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsetup)); 3889566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSmooth %d", (int)i)); 3899566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsolve)); 3909566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASResid %d", (int)i)); 3919566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventresidual)); 3929566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASInterp %d", (int)i)); 3939566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventinterprestrict)); 3940dd27c6cSPeter Brune } else { 3950298fd71SBarry Smith fas->eventsmoothsetup = 0; 3960298fd71SBarry Smith fas->eventsmoothsolve = 0; 3970298fd71SBarry Smith fas->eventresidual = 0; 3980298fd71SBarry Smith fas->eventinterprestrict = 0; 3990dd27c6cSPeter Brune } 4000dd27c6cSPeter Brune } 4010dd27c6cSPeter Brune } 4020dd27c6cSPeter Brune PetscFunctionReturn(0); 4030dd27c6cSPeter Brune } 4040dd27c6cSPeter Brune 405ab8d36c9SPeter Brune /* 406ab8d36c9SPeter Brune Creates the default smoother type. 407ab8d36c9SPeter Brune 40804d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level. 409ab8d36c9SPeter Brune 410ab8d36c9SPeter Brune */ 4119371c9d4SSatish Balay PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) { 412ab8d36c9SPeter Brune SNES_FAS *fas; 413ab8d36c9SPeter Brune const char *optionsprefix; 414ab8d36c9SPeter Brune char tprefix[128]; 415ab8d36c9SPeter Brune SNES nsmooth; 41622d28d08SBarry Smith 417ab8d36c9SPeter Brune PetscFunctionBegin; 418f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 419064a246eSJacob Faibussowitsch PetscValidPointer(smooth, 2); 420ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 4219566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 422ab8d36c9SPeter Brune /* create the default smoother */ 4239566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth)); 424ab8d36c9SPeter Brune if (fas->level == 0) { 4259566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_coarse_", sizeof(tprefix))); 4269566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix)); 4279566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix)); 4289566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNEWTONLS)); 4299566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs)); 430ab8d36c9SPeter Brune } else { 4319566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%d_", (int)fas->level)); 4329566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix)); 4339566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix)); 4349566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNRICHARDSON)); 4359566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs)); 436ab8d36c9SPeter Brune } 4379566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1)); 4389566063dSJacob Faibussowitsch PetscCall(PetscLogObjectParent((PetscObject)snes, (PetscObject)nsmooth)); 4399566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth)); 4409566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)nsmooth, PetscMGLevelId, fas->level)); 441ab8d36c9SPeter Brune *smooth = nsmooth; 442ab8d36c9SPeter Brune PetscFunctionReturn(0); 443ab8d36c9SPeter Brune } 444ab8d36c9SPeter Brune 445ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */ 446ab8d36c9SPeter Brune 447ab8d36c9SPeter Brune /*@ 448ab8d36c9SPeter Brune SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 449ab8d36c9SPeter Brune 450*f6dfbefdSBarry Smith Logically Collective on snes 451ab8d36c9SPeter Brune 452ab8d36c9SPeter Brune Input Parameters: 453*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 454ab8d36c9SPeter Brune . level - the level to set the number of cycles on 455ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 456ab8d36c9SPeter Brune 457ab8d36c9SPeter Brune Level: advanced 458ab8d36c9SPeter Brune 459*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetCycles()` 460ab8d36c9SPeter Brune @*/ 4619371c9d4SSatish Balay PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) { 462f833ba53SLisandro Dalcin SNES_FAS *fas; 46322d28d08SBarry Smith 464ab8d36c9SPeter Brune PetscFunctionBegin; 465f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 466f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 467ab8d36c9SPeter Brune fas->n_cycles = cycles; 4689566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs)); 469ab8d36c9SPeter Brune PetscFunctionReturn(0); 470ab8d36c9SPeter Brune } 471ab8d36c9SPeter Brune 472ab8d36c9SPeter Brune /*@ 473ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 474ab8d36c9SPeter Brune 475*f6dfbefdSBarry Smith Logically Collective on snes 476ab8d36c9SPeter Brune 477*f6dfbefdSBarry Smith Input Parameter: 478*f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 479ab8d36c9SPeter Brune 480*f6dfbefdSBarry Smith Output Parameter: 481ab8d36c9SPeter Brune . smooth - the smoother 482ab8d36c9SPeter Brune 483ab8d36c9SPeter Brune Level: advanced 484ab8d36c9SPeter Brune 485*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()` 486ab8d36c9SPeter Brune @*/ 4879371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) { 488ab8d36c9SPeter Brune SNES_FAS *fas; 4895fd66863SKarl Rupp 490ab8d36c9SPeter Brune PetscFunctionBegin; 491f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 492f833ba53SLisandro Dalcin PetscValidPointer(smooth, 2); 493ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 494ab8d36c9SPeter Brune *smooth = fas->smoothd; 495ab8d36c9SPeter Brune PetscFunctionReturn(0); 496ab8d36c9SPeter Brune } 497ab8d36c9SPeter Brune /*@ 498ab8d36c9SPeter Brune SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 499ab8d36c9SPeter Brune 500*f6dfbefdSBarry Smith Logically Collective on snes 501ab8d36c9SPeter Brune 502*f6dfbefdSBarry Smith Input Parameter: 503*f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 504ab8d36c9SPeter Brune 505*f6dfbefdSBarry Smith Output Parameter: 506ab8d36c9SPeter Brune . smoothu - the smoother 507ab8d36c9SPeter Brune 508*f6dfbefdSBarry Smith Note: 509ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent 510ab8d36c9SPeter Brune default behavior in the process of the solve. 511ab8d36c9SPeter Brune 512ab8d36c9SPeter Brune Level: advanced 513ab8d36c9SPeter Brune 514*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherDown()` 515ab8d36c9SPeter Brune @*/ 5169371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) { 517ab8d36c9SPeter Brune SNES_FAS *fas; 5185fd66863SKarl Rupp 519ab8d36c9SPeter Brune PetscFunctionBegin; 520f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 521f833ba53SLisandro Dalcin PetscValidPointer(smoothu, 2); 522ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 5231aa26658SKarl Rupp if (!fas->smoothu) *smoothu = fas->smoothd; 5241aa26658SKarl Rupp else *smoothu = fas->smoothu; 525ab8d36c9SPeter Brune PetscFunctionReturn(0); 526ab8d36c9SPeter Brune } 527ab8d36c9SPeter Brune 528ab8d36c9SPeter Brune /*@ 529ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 530ab8d36c9SPeter Brune 531*f6dfbefdSBarry Smith Logically Collective on snes 532ab8d36c9SPeter Brune 533*f6dfbefdSBarry Smith Input Parameter: 534*f6dfbefdSBarry Smith . snes - `SNESFAS`, the nonlinear multigrid context 535ab8d36c9SPeter Brune 536*f6dfbefdSBarry Smith Output Parameter: 537ab8d36c9SPeter Brune . smoothd - the smoother 538ab8d36c9SPeter Brune 539ab8d36c9SPeter Brune Level: advanced 540ab8d36c9SPeter Brune 541*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 542ab8d36c9SPeter Brune @*/ 5439371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) { 544ab8d36c9SPeter Brune SNES_FAS *fas; 5455fd66863SKarl Rupp 546ab8d36c9SPeter Brune PetscFunctionBegin; 547f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 548f833ba53SLisandro Dalcin PetscValidPointer(smoothd, 2); 549ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 550ab8d36c9SPeter Brune *smoothd = fas->smoothd; 551ab8d36c9SPeter Brune PetscFunctionReturn(0); 552ab8d36c9SPeter Brune } 553ab8d36c9SPeter Brune 554ab8d36c9SPeter Brune /*@ 555ab8d36c9SPeter Brune SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 556ab8d36c9SPeter Brune 557*f6dfbefdSBarry Smith Logically Collective on snes 558ab8d36c9SPeter Brune 559*f6dfbefdSBarry Smith Input Parameter: 560*f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 561ab8d36c9SPeter Brune 562*f6dfbefdSBarry Smith Output Parameter: 563*f6dfbefdSBarry Smith . correction - the coarse correction solve on this level 564ab8d36c9SPeter Brune 565*f6dfbefdSBarry Smith Note: 5660298fd71SBarry Smith Returns NULL on the coarsest level. 567ab8d36c9SPeter Brune 568ab8d36c9SPeter Brune Level: advanced 569ab8d36c9SPeter Brune 570*f6dfbefdSBarry Smith .seealso: `SNESFAS` `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 571ab8d36c9SPeter Brune @*/ 5729371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) { 573ab8d36c9SPeter Brune SNES_FAS *fas; 5745fd66863SKarl Rupp 575ab8d36c9SPeter Brune PetscFunctionBegin; 576f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 577f833ba53SLisandro Dalcin PetscValidPointer(correction, 2); 578ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 579ab8d36c9SPeter Brune *correction = fas->next; 580ab8d36c9SPeter Brune PetscFunctionReturn(0); 581ab8d36c9SPeter Brune } 582ab8d36c9SPeter Brune 583ab8d36c9SPeter Brune /*@ 584ab8d36c9SPeter Brune SNESFASCycleGetInterpolation - Gets the interpolation on this level 585ab8d36c9SPeter Brune 586*f6dfbefdSBarry Smith Logically Collective on snes 587ab8d36c9SPeter Brune 588*f6dfbefdSBarry Smith Input Parameter: 589*f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 590ab8d36c9SPeter Brune 591*f6dfbefdSBarry Smith Output Parameter: 592ab8d36c9SPeter Brune . mat - the interpolation operator on this level 593ab8d36c9SPeter Brune 594*f6dfbefdSBarry Smith Level: advanced 595ab8d36c9SPeter Brune 596*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 597ab8d36c9SPeter Brune @*/ 5989371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) { 599ab8d36c9SPeter Brune SNES_FAS *fas; 6005fd66863SKarl Rupp 601ab8d36c9SPeter Brune PetscFunctionBegin; 602f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 603f833ba53SLisandro Dalcin PetscValidPointer(mat, 2); 604ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 605ab8d36c9SPeter Brune *mat = fas->interpolate; 606ab8d36c9SPeter Brune PetscFunctionReturn(0); 607ab8d36c9SPeter Brune } 608ab8d36c9SPeter Brune 609ab8d36c9SPeter Brune /*@ 610ab8d36c9SPeter Brune SNESFASCycleGetRestriction - Gets the restriction on this level 611ab8d36c9SPeter Brune 612*f6dfbefdSBarry Smith Logically Collective on snes 613ab8d36c9SPeter Brune 614*f6dfbefdSBarry Smith Input Parameter: 615*f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 616ab8d36c9SPeter Brune 617*f6dfbefdSBarry Smith Output Parameter: 618ab8d36c9SPeter Brune . mat - the restriction operator on this level 619ab8d36c9SPeter Brune 620*f6dfbefdSBarry Smith Level: advanced 621ab8d36c9SPeter Brune 622*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetRestriction()`, `SNESFASCycleGetInterpolation()` 623ab8d36c9SPeter Brune @*/ 6249371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) { 625ab8d36c9SPeter Brune SNES_FAS *fas; 6265fd66863SKarl Rupp 627ab8d36c9SPeter Brune PetscFunctionBegin; 628f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 629f833ba53SLisandro Dalcin PetscValidPointer(mat, 2); 630ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 631ab8d36c9SPeter Brune *mat = fas->restrct; 632ab8d36c9SPeter Brune PetscFunctionReturn(0); 633ab8d36c9SPeter Brune } 634ab8d36c9SPeter Brune 635ab8d36c9SPeter Brune /*@ 636ab8d36c9SPeter Brune SNESFASCycleGetInjection - Gets the injection on this level 637ab8d36c9SPeter Brune 638*f6dfbefdSBarry Smith Logically Collective on snes 639ab8d36c9SPeter Brune 640*f6dfbefdSBarry Smith Input Parameter: 641*f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 642ab8d36c9SPeter Brune 643*f6dfbefdSBarry Smith Output Parameter: 644ab8d36c9SPeter Brune . mat - the restriction operator on this level 645ab8d36c9SPeter Brune 646*f6dfbefdSBarry Smith Level: advanced 647ab8d36c9SPeter Brune 648*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetInjection()`, `SNESFASCycleGetRestriction()` 649ab8d36c9SPeter Brune @*/ 6509371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) { 651ab8d36c9SPeter Brune SNES_FAS *fas; 6525fd66863SKarl Rupp 653ab8d36c9SPeter Brune PetscFunctionBegin; 654f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 655f833ba53SLisandro Dalcin PetscValidPointer(mat, 2); 656ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 657ab8d36c9SPeter Brune *mat = fas->inject; 658ab8d36c9SPeter Brune PetscFunctionReturn(0); 659ab8d36c9SPeter Brune } 660ab8d36c9SPeter Brune 661ab8d36c9SPeter Brune /*@ 662ab8d36c9SPeter Brune SNESFASCycleGetRScale - Gets the injection on this level 663ab8d36c9SPeter Brune 664*f6dfbefdSBarry Smith Logically Collective on snes 665ab8d36c9SPeter Brune 666*f6dfbefdSBarry Smith Input Parameter: 667*f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 668ab8d36c9SPeter Brune 669*f6dfbefdSBarry Smith Output Parameter: 670ab8d36c9SPeter Brune . mat - the restriction operator on this level 671ab8d36c9SPeter Brune 672*f6dfbefdSBarry Smith Level: advanced 673ab8d36c9SPeter Brune 674*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetRestriction()`, `SNESFASGetRScale()` 675ab8d36c9SPeter Brune @*/ 6769371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) { 677ab8d36c9SPeter Brune SNES_FAS *fas; 6785fd66863SKarl Rupp 679ab8d36c9SPeter Brune PetscFunctionBegin; 680f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 681f833ba53SLisandro Dalcin PetscValidPointer(vec, 2); 682ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 683ab8d36c9SPeter Brune *vec = fas->rscale; 684ab8d36c9SPeter Brune PetscFunctionReturn(0); 685ab8d36c9SPeter Brune } 686ab8d36c9SPeter Brune 687ab8d36c9SPeter Brune /*@ 688ab8d36c9SPeter Brune SNESFASCycleIsFine - Determines if a given cycle is the fine level. 689ab8d36c9SPeter Brune 690*f6dfbefdSBarry Smith Logically Collective on snes 691ab8d36c9SPeter Brune 692*f6dfbefdSBarry Smith Input Parameter: 693*f6dfbefdSBarry Smith . snes - the `SNESFAS` `SNES` context 694ab8d36c9SPeter Brune 695*f6dfbefdSBarry Smith Output Parameter: 696ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not 697ab8d36c9SPeter Brune 698ab8d36c9SPeter Brune Level: advanced 699ab8d36c9SPeter Brune 700*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()` 701ab8d36c9SPeter Brune @*/ 7029371c9d4SSatish Balay PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) { 703ab8d36c9SPeter Brune SNES_FAS *fas; 7045fd66863SKarl Rupp 705ab8d36c9SPeter Brune PetscFunctionBegin; 706f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 707534a8f05SLisandro Dalcin PetscValidBoolPointer(flg, 2); 708ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 7091aa26658SKarl Rupp if (fas->level == fas->levels - 1) *flg = PETSC_TRUE; 7101aa26658SKarl Rupp else *flg = PETSC_FALSE; 711ab8d36c9SPeter Brune PetscFunctionReturn(0); 712ab8d36c9SPeter Brune } 713ab8d36c9SPeter Brune 714*f6dfbefdSBarry Smith /* functions called on the finest level that return level-specific information */ 715ab8d36c9SPeter Brune 716ab8d36c9SPeter Brune /*@ 717*f6dfbefdSBarry Smith SNESFASSetInterpolation - Sets the `Mat` to be used to apply the 718ab8d36c9SPeter Brune interpolation from l-1 to the lth level 719ab8d36c9SPeter Brune 720ab8d36c9SPeter Brune Input Parameters: 721*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 722ab8d36c9SPeter Brune . mat - the interpolation operator 723ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 724ab8d36c9SPeter Brune 725ab8d36c9SPeter Brune Level: advanced 726ab8d36c9SPeter Brune 727ab8d36c9SPeter Brune Notes: 728ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction 729ab8d36c9SPeter Brune for the same level. 730ab8d36c9SPeter Brune 731ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 732ab8d36c9SPeter Brune out from the matrix size which one it is. 733ab8d36c9SPeter Brune 734*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`, `SNESFASSetRScale()` 735ab8d36c9SPeter Brune @*/ 7369371c9d4SSatish Balay PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 73722d28d08SBarry Smith SNES_FAS *fas; 738ab8d36c9SPeter Brune SNES levelsnes; 73922d28d08SBarry Smith 740ab8d36c9SPeter Brune PetscFunctionBegin; 741f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 742f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 7439566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 744ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 7459566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 7469566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->interpolate)); 747ab8d36c9SPeter Brune fas->interpolate = mat; 748ab8d36c9SPeter Brune PetscFunctionReturn(0); 749ab8d36c9SPeter Brune } 750ab8d36c9SPeter Brune 751ab8d36c9SPeter Brune /*@ 752ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the 753ab8d36c9SPeter Brune interpolation from l-1 to the lth level 754ab8d36c9SPeter Brune 755ab8d36c9SPeter Brune Input Parameters: 756*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 757ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 758ab8d36c9SPeter Brune 759*f6dfbefdSBarry Smith Output Parameter: 760ab8d36c9SPeter Brune . mat - the interpolation operator 761ab8d36c9SPeter Brune 762ab8d36c9SPeter Brune Level: advanced 763ab8d36c9SPeter Brune 764*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASGetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetRScale()` 765ab8d36c9SPeter Brune @*/ 7669371c9d4SSatish Balay PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) { 76722d28d08SBarry Smith SNES_FAS *fas; 768ab8d36c9SPeter Brune SNES levelsnes; 76922d28d08SBarry Smith 770ab8d36c9SPeter Brune PetscFunctionBegin; 771f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 772f833ba53SLisandro Dalcin PetscValidPointer(mat, 3); 7739566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 774ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 775ab8d36c9SPeter Brune *mat = fas->interpolate; 776ab8d36c9SPeter Brune PetscFunctionReturn(0); 777ab8d36c9SPeter Brune } 778ab8d36c9SPeter Brune 779ab8d36c9SPeter Brune /*@ 780*f6dfbefdSBarry Smith SNESFASSetRestriction - Sets the matrix to be used to restrict the defect 781ab8d36c9SPeter Brune from level l to l-1. 782ab8d36c9SPeter Brune 783ab8d36c9SPeter Brune Input Parameters: 784*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 785ab8d36c9SPeter Brune . mat - the restriction matrix 786ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 787ab8d36c9SPeter Brune 788ab8d36c9SPeter Brune Level: advanced 789ab8d36c9SPeter Brune 790ab8d36c9SPeter Brune Notes: 791ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation 792ab8d36c9SPeter Brune for the same level. 793ab8d36c9SPeter Brune 794ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 795ab8d36c9SPeter Brune out from the matrix size which one it is. 796ab8d36c9SPeter Brune 797ab8d36c9SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 798ab8d36c9SPeter Brune is used. 799ab8d36c9SPeter Brune 800*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetInjection()` 801ab8d36c9SPeter Brune @*/ 8029371c9d4SSatish Balay PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 80322d28d08SBarry Smith SNES_FAS *fas; 804ab8d36c9SPeter Brune SNES levelsnes; 80522d28d08SBarry Smith 806ab8d36c9SPeter Brune PetscFunctionBegin; 807f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 808f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 8099566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 810ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 8119566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 8129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->restrct)); 813ab8d36c9SPeter Brune fas->restrct = mat; 814ab8d36c9SPeter Brune PetscFunctionReturn(0); 815ab8d36c9SPeter Brune } 816ab8d36c9SPeter Brune 817ab8d36c9SPeter Brune /*@ 818ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the 819ab8d36c9SPeter Brune restriction from l to the l-1th level 820ab8d36c9SPeter Brune 821ab8d36c9SPeter Brune Input Parameters: 822*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 823ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 824ab8d36c9SPeter Brune 825*f6dfbefdSBarry Smith Output Parameter: 826ab8d36c9SPeter Brune . mat - the interpolation operator 827ab8d36c9SPeter Brune 828ab8d36c9SPeter Brune Level: advanced 829ab8d36c9SPeter Brune 830*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetRestriction()`, `SNESFASGetInjection()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()` 831ab8d36c9SPeter Brune @*/ 8329371c9d4SSatish Balay PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) { 83322d28d08SBarry Smith SNES_FAS *fas; 834ab8d36c9SPeter Brune SNES levelsnes; 83522d28d08SBarry Smith 836ab8d36c9SPeter Brune PetscFunctionBegin; 837f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 838f833ba53SLisandro Dalcin PetscValidPointer(mat, 3); 8399566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 840ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 841ab8d36c9SPeter Brune *mat = fas->restrct; 842ab8d36c9SPeter Brune PetscFunctionReturn(0); 843ab8d36c9SPeter Brune } 844ab8d36c9SPeter Brune 845ab8d36c9SPeter Brune /*@ 846ab8d36c9SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 847ab8d36c9SPeter Brune from level l to l-1. 848ab8d36c9SPeter Brune 849ab8d36c9SPeter Brune Input Parameters: 850*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 851ab8d36c9SPeter Brune . mat - the restriction matrix 852ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 853ab8d36c9SPeter Brune 854ab8d36c9SPeter Brune Level: advanced 855ab8d36c9SPeter Brune 856*f6dfbefdSBarry Smith Note: 857ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to 858ab8d36c9SPeter Brune project the solution instead. 859ab8d36c9SPeter Brune 860*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetRestriction()` 861ab8d36c9SPeter Brune @*/ 8629371c9d4SSatish Balay PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 86322d28d08SBarry Smith SNES_FAS *fas; 864ab8d36c9SPeter Brune SNES levelsnes; 86522d28d08SBarry Smith 866ab8d36c9SPeter Brune PetscFunctionBegin; 867f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 868f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 8699566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 870ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 8719566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 8729566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->inject)); 8731aa26658SKarl Rupp 874ab8d36c9SPeter Brune fas->inject = mat; 875ab8d36c9SPeter Brune PetscFunctionReturn(0); 876ab8d36c9SPeter Brune } 877ab8d36c9SPeter Brune 878ab8d36c9SPeter Brune /*@ 879ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the 880ab8d36c9SPeter Brune injection from l-1 to the lth level 881ab8d36c9SPeter Brune 882ab8d36c9SPeter Brune Input Parameters: 883*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 884ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 885ab8d36c9SPeter Brune 886*f6dfbefdSBarry Smith Output Parameter: 887ab8d36c9SPeter Brune . mat - the injection operator 888ab8d36c9SPeter Brune 889ab8d36c9SPeter Brune Level: advanced 890ab8d36c9SPeter Brune 891*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()` 892ab8d36c9SPeter Brune @*/ 8939371c9d4SSatish Balay PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) { 89422d28d08SBarry Smith SNES_FAS *fas; 895ab8d36c9SPeter Brune SNES levelsnes; 89622d28d08SBarry Smith 897ab8d36c9SPeter Brune PetscFunctionBegin; 898f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 899f833ba53SLisandro Dalcin PetscValidPointer(mat, 3); 9009566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 901ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 902ab8d36c9SPeter Brune *mat = fas->inject; 903ab8d36c9SPeter Brune PetscFunctionReturn(0); 904ab8d36c9SPeter Brune } 905ab8d36c9SPeter Brune 906ab8d36c9SPeter Brune /*@ 907ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 908ab8d36c9SPeter Brune operator from level l to l-1. 909ab8d36c9SPeter Brune 910ab8d36c9SPeter Brune Input Parameters: 911*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 912ab8d36c9SPeter Brune . rscale - the restriction scaling 913ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 914ab8d36c9SPeter Brune 915ab8d36c9SPeter Brune Level: advanced 916ab8d36c9SPeter Brune 917*f6dfbefdSBarry Smith Note: 918ab8d36c9SPeter Brune This is only used in the case that the injection is not set. 919ab8d36c9SPeter Brune 920*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 921ab8d36c9SPeter Brune @*/ 9229371c9d4SSatish Balay PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 923ab8d36c9SPeter Brune SNES_FAS *fas; 924ab8d36c9SPeter Brune SNES levelsnes; 92522d28d08SBarry Smith 926ab8d36c9SPeter Brune PetscFunctionBegin; 927f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 928f833ba53SLisandro Dalcin if (rscale) PetscValidHeaderSpecific(rscale, VEC_CLASSID, 3); 9299566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 930ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 9319566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rscale)); 9329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->rscale)); 933ab8d36c9SPeter Brune fas->rscale = rscale; 934ab8d36c9SPeter Brune PetscFunctionReturn(0); 935ab8d36c9SPeter Brune } 936ab8d36c9SPeter Brune 937ab8d36c9SPeter Brune /*@ 938ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level. 939ab8d36c9SPeter Brune 940ab8d36c9SPeter Brune Input Parameters: 941*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 942ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 943ab8d36c9SPeter Brune 944*f6dfbefdSBarry Smith Output Parameter: 945ab8d36c9SPeter Brune smooth - the smoother 946ab8d36c9SPeter Brune 947ab8d36c9SPeter Brune Level: advanced 948ab8d36c9SPeter Brune 949*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 950ab8d36c9SPeter Brune @*/ 9519371c9d4SSatish Balay PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) { 952ab8d36c9SPeter Brune SNES_FAS *fas; 953ab8d36c9SPeter Brune SNES levelsnes; 95422d28d08SBarry Smith 955ab8d36c9SPeter Brune PetscFunctionBegin; 956f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 957f833ba53SLisandro Dalcin PetscValidPointer(smooth, 3); 9589566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 959ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 96048a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 961ab8d36c9SPeter Brune *smooth = fas->smoothd; 962ab8d36c9SPeter Brune PetscFunctionReturn(0); 963ab8d36c9SPeter Brune } 964ab8d36c9SPeter Brune 965ab8d36c9SPeter Brune /*@ 966ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level. 967ab8d36c9SPeter Brune 968ab8d36c9SPeter Brune Input Parameters: 969*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 970ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 971ab8d36c9SPeter Brune 972*f6dfbefdSBarry Smith Output Parameter: 973ab8d36c9SPeter Brune smooth - the smoother 974ab8d36c9SPeter Brune 975ab8d36c9SPeter Brune Level: advanced 976ab8d36c9SPeter Brune 977*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 978ab8d36c9SPeter Brune @*/ 9799371c9d4SSatish Balay PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) { 980ab8d36c9SPeter Brune SNES_FAS *fas; 981ab8d36c9SPeter Brune SNES levelsnes; 98222d28d08SBarry Smith 983ab8d36c9SPeter Brune PetscFunctionBegin; 984f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 985f833ba53SLisandro Dalcin PetscValidPointer(smooth, 3); 9869566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 987ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 988ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 98948a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 99048a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu)); 991ab8d36c9SPeter Brune *smooth = fas->smoothd; 992ab8d36c9SPeter Brune PetscFunctionReturn(0); 993ab8d36c9SPeter Brune } 994ab8d36c9SPeter Brune 995ab8d36c9SPeter Brune /*@ 996ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level. 997ab8d36c9SPeter Brune 998ab8d36c9SPeter Brune Input Parameters: 999*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1000ab8d36c9SPeter Brune - level - the level (0 is coarsest) 1001ab8d36c9SPeter Brune 1002*f6dfbefdSBarry Smith Output Parameter: 1003ab8d36c9SPeter Brune smooth - the smoother 1004ab8d36c9SPeter Brune 1005ab8d36c9SPeter Brune Level: advanced 1006ab8d36c9SPeter Brune 1007*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1008ab8d36c9SPeter Brune @*/ 10099371c9d4SSatish Balay PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) { 1010ab8d36c9SPeter Brune SNES_FAS *fas; 1011ab8d36c9SPeter Brune SNES levelsnes; 101222d28d08SBarry Smith 1013ab8d36c9SPeter Brune PetscFunctionBegin; 1014f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1015f833ba53SLisandro Dalcin PetscValidPointer(smooth, 3); 10169566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 1017ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1018ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 101948a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 102048a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu)); 1021ab8d36c9SPeter Brune *smooth = fas->smoothu; 1022ab8d36c9SPeter Brune PetscFunctionReturn(0); 1023ab8d36c9SPeter Brune } 1024d6ad1212SPeter Brune 1025d6ad1212SPeter Brune /*@ 1026d6ad1212SPeter Brune SNESFASGetCoarseSolve - Gets the coarsest solver. 1027d6ad1212SPeter Brune 1028*f6dfbefdSBarry Smith Input Parameter: 1029*f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 1030d6ad1212SPeter Brune 1031*f6dfbefdSBarry Smith Output Parameter: 1032a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver 1033d6ad1212SPeter Brune 1034d6ad1212SPeter Brune Level: advanced 1035d6ad1212SPeter Brune 1036*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1037d6ad1212SPeter Brune @*/ 10389371c9d4SSatish Balay PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse) { 1039d6ad1212SPeter Brune SNES_FAS *fas; 1040d6ad1212SPeter Brune SNES levelsnes; 104122d28d08SBarry Smith 1042d6ad1212SPeter Brune PetscFunctionBegin; 1043f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1044064a246eSJacob Faibussowitsch PetscValidPointer(coarse, 2); 10459566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, 0, &levelsnes)); 1046d6ad1212SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1047d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 104848a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 1049a3a80b83SMatthew G. Knepley *coarse = fas->smoothd; 1050d6ad1212SPeter Brune PetscFunctionReturn(0); 1051d6ad1212SPeter Brune } 1052928e959bSPeter Brune 1053928e959bSPeter Brune /*@ 1054*f6dfbefdSBarry Smith SNESFASFullSetDownSweep - Smooth during the initial downsweep for `SNESFAS` 1055928e959bSPeter Brune 1056*f6dfbefdSBarry Smith Logically Collective on snes 1057928e959bSPeter Brune 1058928e959bSPeter Brune Input Parameters: 1059*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1060928e959bSPeter Brune - swp - whether to downsweep or not 1061928e959bSPeter Brune 1062928e959bSPeter Brune Options Database Key: 1063928e959bSPeter Brune . -snes_fas_full_downsweep - Sets number of pre-smoothing steps 1064928e959bSPeter Brune 1065928e959bSPeter Brune Level: advanced 1066928e959bSPeter Brune 1067*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 1068928e959bSPeter Brune @*/ 10699371c9d4SSatish Balay PetscErrorCode SNESFASFullSetDownSweep(SNES snes, PetscBool swp) { 1070f833ba53SLisandro Dalcin SNES_FAS *fas; 1071928e959bSPeter Brune 1072928e959bSPeter Brune PetscFunctionBegin; 1073f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1074f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 1075928e959bSPeter Brune fas->full_downsweep = swp; 10761baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetDownSweep(fas->next, swp)); 1077928e959bSPeter Brune PetscFunctionReturn(0); 1078928e959bSPeter Brune } 10796dfbafefSToby Isaac 10806dfbafefSToby Isaac /*@ 10816dfbafefSToby Isaac SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 10826dfbafefSToby Isaac 1083*f6dfbefdSBarry Smith Logically Collective on snes 10846dfbafefSToby Isaac 10856dfbafefSToby Isaac Input Parameters: 1086*f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 10876dfbafefSToby Isaac - total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 10886dfbafefSToby Isaac 10896dfbafefSToby Isaac Options Database Key: 109067b8a455SSatish Balay . -snes_fas_full_total - Use total restriction and interpolation on the initial down and up sweeps for the full FAS cycle 10916dfbafefSToby Isaac 10926dfbafefSToby Isaac Level: advanced 10936dfbafefSToby Isaac 1094*f6dfbefdSBarry Smith Note: 1095*f6dfbefdSBarry Smith This option is only significant if the interpolation of a coarse correction (`MatInterpolate()`) is significantly different from total 1096*f6dfbefdSBarry Smith solution interpolation (`DMInterpolateSolution()`). 10976dfbafefSToby Isaac 1098*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()` 10996dfbafefSToby Isaac @*/ 11009371c9d4SSatish Balay PetscErrorCode SNESFASFullSetTotal(SNES snes, PetscBool total) { 11016dfbafefSToby Isaac SNES_FAS *fas; 11026dfbafefSToby Isaac 11036dfbafefSToby Isaac PetscFunctionBegin; 11046dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 11056dfbafefSToby Isaac fas = (SNES_FAS *)snes->data; 11066dfbafefSToby Isaac fas->full_total = total; 11071baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetTotal(fas->next, total)); 11086dfbafefSToby Isaac PetscFunctionReturn(0); 11096dfbafefSToby Isaac } 11106dfbafefSToby Isaac 11116dfbafefSToby Isaac /*@ 11126dfbafefSToby Isaac SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 11136dfbafefSToby Isaac 1114*f6dfbefdSBarry Smith Logically Collective on snes 11156dfbafefSToby Isaac 1116*f6dfbefdSBarry Smith Input Parameter: 1117*f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 11186dfbafefSToby Isaac 11196dfbafefSToby Isaac Output: 11206dfbafefSToby Isaac . total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 11216dfbafefSToby Isaac 11226dfbafefSToby Isaac Level: advanced 11236dfbafefSToby Isaac 1124*f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`, `SNESFullSetTotal()` 11256dfbafefSToby Isaac @*/ 11269371c9d4SSatish Balay PetscErrorCode SNESFASFullGetTotal(SNES snes, PetscBool *total) { 11276dfbafefSToby Isaac SNES_FAS *fas; 11286dfbafefSToby Isaac 11296dfbafefSToby Isaac PetscFunctionBegin; 11306dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 11316dfbafefSToby Isaac fas = (SNES_FAS *)snes->data; 11326dfbafefSToby Isaac *total = fas->full_total; 11336dfbafefSToby Isaac PetscFunctionReturn(0); 11346dfbafefSToby Isaac } 1135