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 10f6dfbefdSBarry Smith - fastype - `SNES_FAS_ADDITIVE`, `SNES_FAS_MULTIPLICATIVE`, `SNES_FAS_FULL`, or `SNES_FAS_KASKADE` 11583a1250SSatish Balay 12583a1250SSatish Balay Level: intermediate 13ab8d36c9SPeter Brune 14f6dfbefdSBarry Smith .seealso: `SNESFAS`, `PCMGSetType()`, `SNESFASGetType()` 15ab8d36c9SPeter Brune @*/ 16d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetType(SNES snes, SNESFASType fastype) 17d71ae5a4SJacob Faibussowitsch { 18f833ba53SLisandro Dalcin SNES_FAS *fas; 1922d28d08SBarry Smith 20ab8d36c9SPeter Brune PetscFunctionBegin; 21f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 22ab8d36c9SPeter Brune PetscValidLogicalCollectiveEnum(snes, fastype, 2); 23f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 24ab8d36c9SPeter Brune fas->fastype = fastype; 251baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetType(fas->next, fastype)); 26ab8d36c9SPeter Brune PetscFunctionReturn(0); 27ab8d36c9SPeter Brune } 28ab8d36c9SPeter Brune 29ab8d36c9SPeter Brune /*@ 30f6dfbefdSBarry Smith SNESFASGetType - Gets the update and correction type used for FAS. 31ab8d36c9SPeter Brune 32ab8d36c9SPeter Brune Logically Collective 33ab8d36c9SPeter Brune 34f6dfbefdSBarry Smith Input Parameter: 35f6dfbefdSBarry Smith . snes - `SNESFAS` context 36ab8d36c9SPeter Brune 37f6dfbefdSBarry Smith Output Parameter: 38f6dfbefdSBarry Smith . fastype - `SNES_FAS_ADDITIVE` or `SNES_FAS_MULTIPLICATIVE` 39ab8d36c9SPeter Brune 40583a1250SSatish Balay Level: intermediate 41583a1250SSatish Balay 42f6dfbefdSBarry Smith .seealso: `SNESFAS`, `PCMGSetType()`, `SNESFASSetType()` 43ab8d36c9SPeter Brune @*/ 44d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetType(SNES snes, SNESFASType *fastype) 45d71ae5a4SJacob Faibussowitsch { 46f833ba53SLisandro Dalcin SNES_FAS *fas; 47ab8d36c9SPeter Brune 48ab8d36c9SPeter Brune PetscFunctionBegin; 49f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 50ab8d36c9SPeter Brune PetscValidPointer(fastype, 2); 51f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 52ab8d36c9SPeter Brune *fastype = fas->fastype; 53ab8d36c9SPeter Brune PetscFunctionReturn(0); 54ab8d36c9SPeter Brune } 55ab8d36c9SPeter Brune 56ab8d36c9SPeter Brune /*@C 57f6dfbefdSBarry Smith SNESFASSetLevels - Sets the number of levels to use with `SNESFAS`. 58ab8d36c9SPeter Brune Must be called before any other FAS routine. 59ab8d36c9SPeter Brune 60ab8d36c9SPeter Brune Input Parameters: 61ab8d36c9SPeter Brune + snes - the snes context 62ab8d36c9SPeter Brune . levels - the number of levels 63ab8d36c9SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 642bf68e3eSBarry Smith problems on smaller sets of processors. 65ab8d36c9SPeter Brune 66ab8d36c9SPeter Brune Level: intermediate 67ab8d36c9SPeter Brune 68f6dfbefdSBarry Smith Note: 69ab8d36c9SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 70ab8d36c9SPeter Brune for setting the level options rather than the -fas_coarse prefix. 71ab8d36c9SPeter Brune 72f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetLevels()` 73ab8d36c9SPeter Brune @*/ 74d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms) 75d71ae5a4SJacob Faibussowitsch { 76ab8d36c9SPeter Brune PetscInt i; 77ab8d36c9SPeter Brune const char *optionsprefix; 78ab8d36c9SPeter Brune char tprefix[128]; 79f833ba53SLisandro Dalcin SNES_FAS *fas; 80ab8d36c9SPeter Brune SNES prevsnes; 81ab8d36c9SPeter Brune MPI_Comm comm; 8222d28d08SBarry Smith 83ab8d36c9SPeter Brune PetscFunctionBegin; 84f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 85f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 869566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)snes, &comm)); 87ab8d36c9SPeter Brune if (levels == fas->levels) { 8822d28d08SBarry Smith if (!comms) PetscFunctionReturn(0); 89ab8d36c9SPeter Brune } 90ab8d36c9SPeter Brune /* user has changed the number of levels; reset */ 91dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, reset); 92ab8d36c9SPeter Brune /* destroy any coarser levels if necessary */ 939566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&fas->next)); 940298fd71SBarry Smith fas->next = NULL; 950298fd71SBarry Smith fas->previous = NULL; 96ab8d36c9SPeter Brune prevsnes = snes; 97ab8d36c9SPeter Brune /* setup the finest level */ 989566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 999566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)snes, PetscMGLevelId, levels - 1)); 100ab8d36c9SPeter Brune for (i = levels - 1; i >= 0; i--) { 101ab8d36c9SPeter Brune if (comms) comm = comms[i]; 102ab8d36c9SPeter Brune fas->level = i; 103ab8d36c9SPeter Brune fas->levels = levels; 104ab8d36c9SPeter Brune fas->fine = snes; 1050298fd71SBarry Smith fas->next = NULL; 106ab8d36c9SPeter Brune if (i > 0) { 1079566063dSJacob Faibussowitsch PetscCall(SNESCreate(comm, &fas->next)); 1089566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 1099566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%d_cycle_", (int)fas->level)); 1109566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->next, optionsprefix)); 1119566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->next, tprefix)); 1129566063dSJacob Faibussowitsch PetscCall(SNESSetType(fas->next, SNESFAS)); 1139566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs)); 1149566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i)); 1159566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)fas->next, PetscMGLevelId, i - 1)); 1161aa26658SKarl Rupp 117ab8d36c9SPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 1181aa26658SKarl Rupp 119ab8d36c9SPeter Brune prevsnes = fas->next; 120ab8d36c9SPeter Brune fas = (SNES_FAS *)prevsnes->data; 121ab8d36c9SPeter Brune } 122ab8d36c9SPeter Brune } 123ab8d36c9SPeter Brune PetscFunctionReturn(0); 124ab8d36c9SPeter Brune } 125ab8d36c9SPeter Brune 126ab8d36c9SPeter Brune /*@ 127f6dfbefdSBarry Smith SNESFASGetLevels - Gets the number of levels in a `SNESFAS`, including fine and coarse grids 128ab8d36c9SPeter Brune 129ab8d36c9SPeter Brune Input Parameter: 130f6dfbefdSBarry Smith . snes - the `SNES` nonlinear solver context 131ab8d36c9SPeter Brune 132ab8d36c9SPeter Brune Output parameter: 133ab8d36c9SPeter Brune . levels - the number of levels 134ab8d36c9SPeter Brune 135ab8d36c9SPeter Brune Level: advanced 136ab8d36c9SPeter Brune 137f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()`, `PCMGGetLevels()` 138ab8d36c9SPeter Brune @*/ 139d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels) 140d71ae5a4SJacob Faibussowitsch { 141f833ba53SLisandro Dalcin SNES_FAS *fas; 1425fd66863SKarl Rupp 143ab8d36c9SPeter Brune PetscFunctionBegin; 144f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 145064a246eSJacob Faibussowitsch PetscValidIntPointer(levels, 2); 146f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 147ab8d36c9SPeter Brune *levels = fas->levels; 148ab8d36c9SPeter Brune PetscFunctionReturn(0); 149ab8d36c9SPeter Brune } 150ab8d36c9SPeter Brune 151ab8d36c9SPeter Brune /*@ 152f6dfbefdSBarry Smith SNESFASGetCycleSNES - Gets the `SNES` corresponding to a particular 153f6dfbefdSBarry Smith level of the `SNESFAS` hierarchy. 154ab8d36c9SPeter Brune 155ab8d36c9SPeter Brune Input Parameters: 156f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context 157f6dfbefdSBarry Smith - level - the level to get 158f6dfbefdSBarry Smith 159f6dfbefdSBarry Smith Output Parameter: 160f6dfbefdSBarry Smith . lsnes - the `SNES` for the requested level 161ab8d36c9SPeter Brune 162ab8d36c9SPeter Brune Level: advanced 163ab8d36c9SPeter Brune 164f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()`, `SNESFASGetLevels()` 165ab8d36c9SPeter Brune @*/ 166d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetCycleSNES(SNES snes, PetscInt level, SNES *lsnes) 167d71ae5a4SJacob Faibussowitsch { 168f833ba53SLisandro Dalcin SNES_FAS *fas; 169ab8d36c9SPeter Brune PetscInt i; 170ab8d36c9SPeter Brune 171ab8d36c9SPeter Brune PetscFunctionBegin; 172f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 173f833ba53SLisandro Dalcin PetscValidPointer(lsnes, 3); 174f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 17563a3b9bcSJacob 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); 1760b121fc5SBarry Smith PetscCheck(fas->level == fas->levels - 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetCycleSNES may only be called on the finest-level SNES"); 177ab8d36c9SPeter Brune 178ab8d36c9SPeter Brune *lsnes = snes; 179ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) { 180ab8d36c9SPeter Brune *lsnes = fas->next; 181ab8d36c9SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 182ab8d36c9SPeter Brune } 18308401ef6SPierre Jolivet PetscCheck(fas->level == level, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "SNESFAS level hierarchy corrupt"); 184ab8d36c9SPeter Brune PetscFunctionReturn(0); 185ab8d36c9SPeter Brune } 186ab8d36c9SPeter Brune 187ab8d36c9SPeter Brune /*@ 188ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 189ab8d36c9SPeter Brune use on all levels. 190ab8d36c9SPeter Brune 191*c3339decSBarry Smith Logically Collective 192ab8d36c9SPeter Brune 193ab8d36c9SPeter Brune Input Parameters: 194f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context 195ab8d36c9SPeter Brune - n - the number of smoothing steps 196ab8d36c9SPeter Brune 197ab8d36c9SPeter Brune Options Database Key: 198ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 199ab8d36c9SPeter Brune 200ab8d36c9SPeter Brune Level: advanced 201ab8d36c9SPeter Brune 202f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothDown()` 203ab8d36c9SPeter Brune @*/ 204d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) 205d71ae5a4SJacob Faibussowitsch { 206f833ba53SLisandro Dalcin SNES_FAS *fas; 20722d28d08SBarry Smith 208ab8d36c9SPeter Brune PetscFunctionBegin; 209f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 210f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 211ab8d36c9SPeter Brune fas->max_up_it = n; 21248a46eb9SPierre Jolivet if (!fas->smoothu && fas->level != 0) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu)); 2131baa6e33SBarry Smith if (fas->smoothu) PetscCall(SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs)); 2141baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothUp(fas->next, n)); 215ab8d36c9SPeter Brune PetscFunctionReturn(0); 216ab8d36c9SPeter Brune } 217ab8d36c9SPeter Brune 218ab8d36c9SPeter Brune /*@ 219ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 220ab8d36c9SPeter Brune use on all levels. 221ab8d36c9SPeter Brune 222*c3339decSBarry Smith Logically Collective 223ab8d36c9SPeter Brune 224ab8d36c9SPeter Brune Input Parameters: 225f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 226ab8d36c9SPeter Brune - n - the number of smoothing steps 227ab8d36c9SPeter Brune 228ab8d36c9SPeter Brune Options Database Key: 229ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 230ab8d36c9SPeter Brune 231ab8d36c9SPeter Brune Level: advanced 232ab8d36c9SPeter Brune 233f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 234ab8d36c9SPeter Brune @*/ 235d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) 236d71ae5a4SJacob Faibussowitsch { 237f833ba53SLisandro Dalcin SNES_FAS *fas; 23822d28d08SBarry Smith 239ab8d36c9SPeter Brune PetscFunctionBegin; 240f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 241f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 24248a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd)); 2439566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs)); 2441aa26658SKarl Rupp 245ab8d36c9SPeter Brune fas->max_down_it = n; 2461baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothDown(fas->next, n)); 247ab8d36c9SPeter Brune PetscFunctionReturn(0); 248ab8d36c9SPeter Brune } 249ab8d36c9SPeter Brune 25087f44e3fSPeter Brune /*@ 251f6dfbefdSBarry Smith SNESFASSetContinuation - Sets the `SNESFAS` cycle to default to exact Newton solves on the upsweep 25287f44e3fSPeter Brune 253*c3339decSBarry Smith Logically Collective 25487f44e3fSPeter Brune 25587f44e3fSPeter Brune Input Parameters: 256f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 25787f44e3fSPeter Brune - n - the number of smoothing steps 25887f44e3fSPeter Brune 25987f44e3fSPeter Brune Options Database Key: 26087f44e3fSPeter Brune . -snes_fas_continuation - sets continuation to true 26187f44e3fSPeter Brune 26287f44e3fSPeter Brune Level: advanced 26387f44e3fSPeter Brune 264f6dfbefdSBarry Smith Note: 26595452b02SPatrick Sanan This sets the prefix on the upsweep smoothers to -fas_continuation 26687f44e3fSPeter Brune 267f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 26887f44e3fSPeter Brune @*/ 269d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetContinuation(SNES snes, PetscBool continuation) 270d71ae5a4SJacob Faibussowitsch { 27187f44e3fSPeter Brune const char *optionsprefix; 27287f44e3fSPeter Brune char tprefix[128]; 273f833ba53SLisandro Dalcin SNES_FAS *fas; 27487f44e3fSPeter Brune 27587f44e3fSPeter Brune PetscFunctionBegin; 276f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 277f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 2789566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 27948a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu)); 2809566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_levels_continuation_", sizeof(tprefix))); 2819566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(fas->smoothu, optionsprefix)); 2829566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->smoothu, tprefix)); 2839566063dSJacob Faibussowitsch PetscCall(SNESSetType(fas->smoothu, SNESNEWTONLS)); 2849566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothu, fas->fine->abstol, fas->fine->rtol, fas->fine->stol, 50, 100)); 28587f44e3fSPeter Brune fas->continuation = continuation; 2861baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetContinuation(fas->next, continuation)); 28787f44e3fSPeter Brune PetscFunctionReturn(0); 28887f44e3fSPeter Brune } 28987f44e3fSPeter Brune 290ab8d36c9SPeter Brune /*@ 291f6dfbefdSBarry Smith SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use `SNESFASSetCyclesOnLevel()` for more 292ab8d36c9SPeter Brune complicated cycling. 293ab8d36c9SPeter Brune 294*c3339decSBarry Smith Logically Collective 295ab8d36c9SPeter Brune 296ab8d36c9SPeter Brune Input Parameters: 297f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 298ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 299ab8d36c9SPeter Brune 300ab8d36c9SPeter Brune Options Database Key: 30167b8a455SSatish Balay . -snes_fas_cycles <1,2> - 1 for V-cycle, 2 for W-cycle 302ab8d36c9SPeter Brune 303ab8d36c9SPeter Brune Level: advanced 304ab8d36c9SPeter Brune 305f6dfbefdSBarry Smith .seealso: `SNES`, `SNESFAS`, `SNESFASSetCyclesOnLevel()` 306ab8d36c9SPeter Brune @*/ 307d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) 308d71ae5a4SJacob Faibussowitsch { 309f833ba53SLisandro Dalcin SNES_FAS *fas; 310ab8d36c9SPeter Brune PetscBool isFine; 31122d28d08SBarry Smith 312ab8d36c9SPeter Brune PetscFunctionBegin; 313f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3149566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 315f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 316ab8d36c9SPeter Brune fas->n_cycles = cycles; 31748a46eb9SPierre Jolivet if (!isFine) PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs)); 3181baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetCycles(fas->next, cycles)); 319ab8d36c9SPeter Brune PetscFunctionReturn(0); 320ab8d36c9SPeter Brune } 321ab8d36c9SPeter Brune 322c8c899caSPeter Brune /*@ 323c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring 324c8c899caSPeter Brune 325*c3339decSBarry Smith Logically Collective 326c8c899caSPeter Brune 327c8c899caSPeter Brune Input Parameters: 328f6dfbefdSBarry Smith + snes - the `SNESFAS` context 329d142ab34SLawrence Mitchell . vf - viewer and format structure (may be NULL if flg is FALSE) 330c8c899caSPeter Brune - flg - monitor or not 331c8c899caSPeter Brune 332c8c899caSPeter Brune Level: advanced 333c8c899caSPeter Brune 334f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESSetMonitor()`, `SNESFASSetCyclesOnLevel()` 335c8c899caSPeter Brune @*/ 336d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg) 337d71ae5a4SJacob Faibussowitsch { 338f833ba53SLisandro Dalcin SNES_FAS *fas; 339c8c899caSPeter Brune PetscBool isFine; 340f833ba53SLisandro Dalcin PetscInt i, levels; 341c8c899caSPeter Brune SNES levelsnes; 34222d28d08SBarry Smith 343c8c899caSPeter Brune PetscFunctionBegin; 344f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3459566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 346f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 347f833ba53SLisandro Dalcin levels = fas->levels; 348c8c899caSPeter Brune if (isFine) { 349c8c899caSPeter Brune for (i = 0; i < levels; i++) { 3509566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 351c8c899caSPeter Brune fas = (SNES_FAS *)levelsnes->data; 352c8c899caSPeter Brune if (flg) { 353c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */ 3549566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes)); 355d142ab34SLawrence Mitchell /* Only register destroy on finest level */ 3569566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(levelsnes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))SNESMonitorDefault, vf, (!i ? (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy : NULL))); 3571aa26658SKarl Rupp } else if (i != fas->levels - 1) { 358c8c899caSPeter Brune /* unset the monitors on the coarse levels */ 3599566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes)); 360c8c899caSPeter Brune } 361c8c899caSPeter Brune } 362c8c899caSPeter Brune } 363c8c899caSPeter Brune PetscFunctionReturn(0); 364c8c899caSPeter Brune } 365c8c899caSPeter Brune 3660dd27c6cSPeter Brune /*@ 367f6dfbefdSBarry Smith SNESFASSetLog - Sets or unsets time logging for various `SNESFAS` stages on all levels 3680dd27c6cSPeter Brune 369*c3339decSBarry Smith Logically Collective 3700dd27c6cSPeter Brune 3710dd27c6cSPeter Brune Input Parameters: 372f6dfbefdSBarry Smith + snes - the `SNESFAS` context 3730dd27c6cSPeter Brune - flg - monitor or not 3740dd27c6cSPeter Brune 3750dd27c6cSPeter Brune Level: advanced 3760dd27c6cSPeter Brune 377f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetMonitor()` 3780dd27c6cSPeter Brune @*/ 379d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) 380d71ae5a4SJacob Faibussowitsch { 381f833ba53SLisandro Dalcin SNES_FAS *fas; 3820dd27c6cSPeter Brune PetscBool isFine; 383f833ba53SLisandro Dalcin PetscInt i, levels; 3840dd27c6cSPeter Brune SNES levelsnes; 3850dd27c6cSPeter Brune char eventname[128]; 3860dd27c6cSPeter Brune 3870dd27c6cSPeter Brune PetscFunctionBegin; 388f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3899566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 390f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 391f833ba53SLisandro Dalcin levels = fas->levels; 3920dd27c6cSPeter Brune if (isFine) { 3930dd27c6cSPeter Brune for (i = 0; i < levels; i++) { 3949566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 3950dd27c6cSPeter Brune fas = (SNES_FAS *)levelsnes->data; 3960dd27c6cSPeter Brune if (flg) { 3979566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSetup %d", (int)i)); 3989566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsetup)); 3999566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSmooth %d", (int)i)); 4009566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsolve)); 4019566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASResid %d", (int)i)); 4029566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventresidual)); 4039566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASInterp %d", (int)i)); 4049566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventinterprestrict)); 4050dd27c6cSPeter Brune } else { 4060298fd71SBarry Smith fas->eventsmoothsetup = 0; 4070298fd71SBarry Smith fas->eventsmoothsolve = 0; 4080298fd71SBarry Smith fas->eventresidual = 0; 4090298fd71SBarry Smith fas->eventinterprestrict = 0; 4100dd27c6cSPeter Brune } 4110dd27c6cSPeter Brune } 4120dd27c6cSPeter Brune } 4130dd27c6cSPeter Brune PetscFunctionReturn(0); 4140dd27c6cSPeter Brune } 4150dd27c6cSPeter Brune 416ab8d36c9SPeter Brune /* 417ab8d36c9SPeter Brune Creates the default smoother type. 418ab8d36c9SPeter Brune 41904d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level. 420ab8d36c9SPeter Brune 421ab8d36c9SPeter Brune */ 422d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) 423d71ae5a4SJacob Faibussowitsch { 424ab8d36c9SPeter Brune SNES_FAS *fas; 425ab8d36c9SPeter Brune const char *optionsprefix; 426ab8d36c9SPeter Brune char tprefix[128]; 427ab8d36c9SPeter Brune SNES nsmooth; 42822d28d08SBarry Smith 429ab8d36c9SPeter Brune PetscFunctionBegin; 430f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 431064a246eSJacob Faibussowitsch PetscValidPointer(smooth, 2); 432ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 4339566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 434ab8d36c9SPeter Brune /* create the default smoother */ 4359566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth)); 436ab8d36c9SPeter Brune if (fas->level == 0) { 4379566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_coarse_", sizeof(tprefix))); 4389566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix)); 4399566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix)); 4409566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNEWTONLS)); 4419566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs)); 442ab8d36c9SPeter Brune } else { 4439566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%d_", (int)fas->level)); 4449566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix)); 4459566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix)); 4469566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNRICHARDSON)); 4479566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs)); 448ab8d36c9SPeter Brune } 4499566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1)); 4509566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth)); 4519566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)nsmooth, PetscMGLevelId, fas->level)); 452ab8d36c9SPeter Brune *smooth = nsmooth; 453ab8d36c9SPeter Brune PetscFunctionReturn(0); 454ab8d36c9SPeter Brune } 455ab8d36c9SPeter Brune 456ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */ 457ab8d36c9SPeter Brune 458ab8d36c9SPeter Brune /*@ 459ab8d36c9SPeter Brune SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 460ab8d36c9SPeter Brune 461*c3339decSBarry Smith Logically Collective 462ab8d36c9SPeter Brune 463ab8d36c9SPeter Brune Input Parameters: 464f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 465ab8d36c9SPeter Brune . level - the level to set the number of cycles on 466ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 467ab8d36c9SPeter Brune 468ab8d36c9SPeter Brune Level: advanced 469ab8d36c9SPeter Brune 470f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetCycles()` 471ab8d36c9SPeter Brune @*/ 472d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) 473d71ae5a4SJacob Faibussowitsch { 474f833ba53SLisandro Dalcin SNES_FAS *fas; 47522d28d08SBarry Smith 476ab8d36c9SPeter Brune PetscFunctionBegin; 477f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 478f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 479ab8d36c9SPeter Brune fas->n_cycles = cycles; 4809566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs)); 481ab8d36c9SPeter Brune PetscFunctionReturn(0); 482ab8d36c9SPeter Brune } 483ab8d36c9SPeter Brune 484ab8d36c9SPeter Brune /*@ 485ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 486ab8d36c9SPeter Brune 487*c3339decSBarry Smith Logically Collective 488ab8d36c9SPeter Brune 489f6dfbefdSBarry Smith Input Parameter: 490f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 491ab8d36c9SPeter Brune 492f6dfbefdSBarry Smith Output Parameter: 493ab8d36c9SPeter Brune . smooth - the smoother 494ab8d36c9SPeter Brune 495ab8d36c9SPeter Brune Level: advanced 496ab8d36c9SPeter Brune 497f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()` 498ab8d36c9SPeter Brune @*/ 499d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 500d71ae5a4SJacob Faibussowitsch { 501ab8d36c9SPeter Brune SNES_FAS *fas; 5025fd66863SKarl Rupp 503ab8d36c9SPeter Brune PetscFunctionBegin; 504f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 505f833ba53SLisandro Dalcin PetscValidPointer(smooth, 2); 506ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 507ab8d36c9SPeter Brune *smooth = fas->smoothd; 508ab8d36c9SPeter Brune PetscFunctionReturn(0); 509ab8d36c9SPeter Brune } 510ab8d36c9SPeter Brune /*@ 511ab8d36c9SPeter Brune SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 512ab8d36c9SPeter Brune 513*c3339decSBarry Smith Logically Collective 514ab8d36c9SPeter Brune 515f6dfbefdSBarry Smith Input Parameter: 516f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 517ab8d36c9SPeter Brune 518f6dfbefdSBarry Smith Output Parameter: 519ab8d36c9SPeter Brune . smoothu - the smoother 520ab8d36c9SPeter Brune 521f6dfbefdSBarry Smith Note: 522ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent 523ab8d36c9SPeter Brune default behavior in the process of the solve. 524ab8d36c9SPeter Brune 525ab8d36c9SPeter Brune Level: advanced 526ab8d36c9SPeter Brune 527f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherDown()` 528ab8d36c9SPeter Brune @*/ 529d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 530d71ae5a4SJacob Faibussowitsch { 531ab8d36c9SPeter Brune SNES_FAS *fas; 5325fd66863SKarl Rupp 533ab8d36c9SPeter Brune PetscFunctionBegin; 534f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 535f833ba53SLisandro Dalcin PetscValidPointer(smoothu, 2); 536ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 5371aa26658SKarl Rupp if (!fas->smoothu) *smoothu = fas->smoothd; 5381aa26658SKarl Rupp else *smoothu = fas->smoothu; 539ab8d36c9SPeter Brune PetscFunctionReturn(0); 540ab8d36c9SPeter Brune } 541ab8d36c9SPeter Brune 542ab8d36c9SPeter Brune /*@ 543ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 544ab8d36c9SPeter Brune 545*c3339decSBarry Smith Logically Collective 546ab8d36c9SPeter Brune 547f6dfbefdSBarry Smith Input Parameter: 548f6dfbefdSBarry Smith . snes - `SNESFAS`, the nonlinear multigrid context 549ab8d36c9SPeter Brune 550f6dfbefdSBarry Smith Output Parameter: 551ab8d36c9SPeter Brune . smoothd - the smoother 552ab8d36c9SPeter Brune 553ab8d36c9SPeter Brune Level: advanced 554ab8d36c9SPeter Brune 555f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 556ab8d36c9SPeter Brune @*/ 557d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 558d71ae5a4SJacob Faibussowitsch { 559ab8d36c9SPeter Brune SNES_FAS *fas; 5605fd66863SKarl Rupp 561ab8d36c9SPeter Brune PetscFunctionBegin; 562f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 563f833ba53SLisandro Dalcin PetscValidPointer(smoothd, 2); 564ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 565ab8d36c9SPeter Brune *smoothd = fas->smoothd; 566ab8d36c9SPeter Brune PetscFunctionReturn(0); 567ab8d36c9SPeter Brune } 568ab8d36c9SPeter Brune 569ab8d36c9SPeter Brune /*@ 570ab8d36c9SPeter Brune SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 571ab8d36c9SPeter Brune 572*c3339decSBarry Smith Logically Collective 573ab8d36c9SPeter Brune 574f6dfbefdSBarry Smith Input Parameter: 575f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 576ab8d36c9SPeter Brune 577f6dfbefdSBarry Smith Output Parameter: 578f6dfbefdSBarry Smith . correction - the coarse correction solve on this level 579ab8d36c9SPeter Brune 580f6dfbefdSBarry Smith Note: 5810298fd71SBarry Smith Returns NULL on the coarsest level. 582ab8d36c9SPeter Brune 583ab8d36c9SPeter Brune Level: advanced 584ab8d36c9SPeter Brune 585f6dfbefdSBarry Smith .seealso: `SNESFAS` `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 586ab8d36c9SPeter Brune @*/ 587d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 588d71ae5a4SJacob Faibussowitsch { 589ab8d36c9SPeter Brune SNES_FAS *fas; 5905fd66863SKarl Rupp 591ab8d36c9SPeter Brune PetscFunctionBegin; 592f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 593f833ba53SLisandro Dalcin PetscValidPointer(correction, 2); 594ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 595ab8d36c9SPeter Brune *correction = fas->next; 596ab8d36c9SPeter Brune PetscFunctionReturn(0); 597ab8d36c9SPeter Brune } 598ab8d36c9SPeter Brune 599ab8d36c9SPeter Brune /*@ 600ab8d36c9SPeter Brune SNESFASCycleGetInterpolation - Gets the interpolation on this level 601ab8d36c9SPeter Brune 602*c3339decSBarry Smith Logically Collective 603ab8d36c9SPeter Brune 604f6dfbefdSBarry Smith Input Parameter: 605f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 606ab8d36c9SPeter Brune 607f6dfbefdSBarry Smith Output Parameter: 608ab8d36c9SPeter Brune . mat - the interpolation operator on this level 609ab8d36c9SPeter Brune 610f6dfbefdSBarry Smith Level: advanced 611ab8d36c9SPeter Brune 612f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 613ab8d36c9SPeter Brune @*/ 614d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 615d71ae5a4SJacob Faibussowitsch { 616ab8d36c9SPeter Brune SNES_FAS *fas; 6175fd66863SKarl Rupp 618ab8d36c9SPeter Brune PetscFunctionBegin; 619f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 620f833ba53SLisandro Dalcin PetscValidPointer(mat, 2); 621ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 622ab8d36c9SPeter Brune *mat = fas->interpolate; 623ab8d36c9SPeter Brune PetscFunctionReturn(0); 624ab8d36c9SPeter Brune } 625ab8d36c9SPeter Brune 626ab8d36c9SPeter Brune /*@ 627ab8d36c9SPeter Brune SNESFASCycleGetRestriction - Gets the restriction on this level 628ab8d36c9SPeter Brune 629*c3339decSBarry Smith Logically Collective 630ab8d36c9SPeter Brune 631f6dfbefdSBarry Smith Input Parameter: 632f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 633ab8d36c9SPeter Brune 634f6dfbefdSBarry Smith Output Parameter: 635ab8d36c9SPeter Brune . mat - the restriction operator on this level 636ab8d36c9SPeter Brune 637f6dfbefdSBarry Smith Level: advanced 638ab8d36c9SPeter Brune 639f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetRestriction()`, `SNESFASCycleGetInterpolation()` 640ab8d36c9SPeter Brune @*/ 641d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 642d71ae5a4SJacob Faibussowitsch { 643ab8d36c9SPeter Brune SNES_FAS *fas; 6445fd66863SKarl Rupp 645ab8d36c9SPeter Brune PetscFunctionBegin; 646f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 647f833ba53SLisandro Dalcin PetscValidPointer(mat, 2); 648ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 649ab8d36c9SPeter Brune *mat = fas->restrct; 650ab8d36c9SPeter Brune PetscFunctionReturn(0); 651ab8d36c9SPeter Brune } 652ab8d36c9SPeter Brune 653ab8d36c9SPeter Brune /*@ 654ab8d36c9SPeter Brune SNESFASCycleGetInjection - Gets the injection on this level 655ab8d36c9SPeter Brune 656*c3339decSBarry Smith Logically Collective 657ab8d36c9SPeter Brune 658f6dfbefdSBarry Smith Input Parameter: 659f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 660ab8d36c9SPeter Brune 661f6dfbefdSBarry Smith Output Parameter: 662ab8d36c9SPeter Brune . mat - the restriction operator on this level 663ab8d36c9SPeter Brune 664f6dfbefdSBarry Smith Level: advanced 665ab8d36c9SPeter Brune 666f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetInjection()`, `SNESFASCycleGetRestriction()` 667ab8d36c9SPeter Brune @*/ 668d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 669d71ae5a4SJacob Faibussowitsch { 670ab8d36c9SPeter Brune SNES_FAS *fas; 6715fd66863SKarl Rupp 672ab8d36c9SPeter Brune PetscFunctionBegin; 673f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 674f833ba53SLisandro Dalcin PetscValidPointer(mat, 2); 675ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 676ab8d36c9SPeter Brune *mat = fas->inject; 677ab8d36c9SPeter Brune PetscFunctionReturn(0); 678ab8d36c9SPeter Brune } 679ab8d36c9SPeter Brune 680ab8d36c9SPeter Brune /*@ 681ab8d36c9SPeter Brune SNESFASCycleGetRScale - Gets the injection on this level 682ab8d36c9SPeter Brune 683*c3339decSBarry Smith Logically Collective 684ab8d36c9SPeter Brune 685f6dfbefdSBarry Smith Input Parameter: 686f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 687ab8d36c9SPeter Brune 688f6dfbefdSBarry Smith Output Parameter: 689ab8d36c9SPeter Brune . mat - the restriction operator on this level 690ab8d36c9SPeter Brune 691f6dfbefdSBarry Smith Level: advanced 692ab8d36c9SPeter Brune 693f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetRestriction()`, `SNESFASGetRScale()` 694ab8d36c9SPeter Brune @*/ 695d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 696d71ae5a4SJacob Faibussowitsch { 697ab8d36c9SPeter Brune SNES_FAS *fas; 6985fd66863SKarl Rupp 699ab8d36c9SPeter Brune PetscFunctionBegin; 700f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 701f833ba53SLisandro Dalcin PetscValidPointer(vec, 2); 702ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 703ab8d36c9SPeter Brune *vec = fas->rscale; 704ab8d36c9SPeter Brune PetscFunctionReturn(0); 705ab8d36c9SPeter Brune } 706ab8d36c9SPeter Brune 707ab8d36c9SPeter Brune /*@ 708ab8d36c9SPeter Brune SNESFASCycleIsFine - Determines if a given cycle is the fine level. 709ab8d36c9SPeter Brune 710*c3339decSBarry Smith Logically Collective 711ab8d36c9SPeter Brune 712f6dfbefdSBarry Smith Input Parameter: 713f6dfbefdSBarry Smith . snes - the `SNESFAS` `SNES` context 714ab8d36c9SPeter Brune 715f6dfbefdSBarry Smith Output Parameter: 716ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not 717ab8d36c9SPeter Brune 718ab8d36c9SPeter Brune Level: advanced 719ab8d36c9SPeter Brune 720f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()` 721ab8d36c9SPeter Brune @*/ 722d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 723d71ae5a4SJacob Faibussowitsch { 724ab8d36c9SPeter Brune SNES_FAS *fas; 7255fd66863SKarl Rupp 726ab8d36c9SPeter Brune PetscFunctionBegin; 727f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 728534a8f05SLisandro Dalcin PetscValidBoolPointer(flg, 2); 729ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 7301aa26658SKarl Rupp if (fas->level == fas->levels - 1) *flg = PETSC_TRUE; 7311aa26658SKarl Rupp else *flg = PETSC_FALSE; 732ab8d36c9SPeter Brune PetscFunctionReturn(0); 733ab8d36c9SPeter Brune } 734ab8d36c9SPeter Brune 735f6dfbefdSBarry Smith /* functions called on the finest level that return level-specific information */ 736ab8d36c9SPeter Brune 737ab8d36c9SPeter Brune /*@ 738f6dfbefdSBarry Smith SNESFASSetInterpolation - Sets the `Mat` to be used to apply the 739ab8d36c9SPeter Brune interpolation from l-1 to the lth level 740ab8d36c9SPeter Brune 741ab8d36c9SPeter Brune Input Parameters: 742f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 743ab8d36c9SPeter Brune . mat - the interpolation operator 744ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 745ab8d36c9SPeter Brune 746ab8d36c9SPeter Brune Level: advanced 747ab8d36c9SPeter Brune 748ab8d36c9SPeter Brune Notes: 749ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction 750ab8d36c9SPeter Brune for the same level. 751ab8d36c9SPeter Brune 752ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 753ab8d36c9SPeter Brune out from the matrix size which one it is. 754ab8d36c9SPeter Brune 755f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`, `SNESFASSetRScale()` 756ab8d36c9SPeter Brune @*/ 757d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) 758d71ae5a4SJacob Faibussowitsch { 75922d28d08SBarry Smith SNES_FAS *fas; 760ab8d36c9SPeter Brune SNES levelsnes; 76122d28d08SBarry Smith 762ab8d36c9SPeter Brune PetscFunctionBegin; 763f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 764f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 7659566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 766ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 7679566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 7689566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->interpolate)); 769ab8d36c9SPeter Brune fas->interpolate = mat; 770ab8d36c9SPeter Brune PetscFunctionReturn(0); 771ab8d36c9SPeter Brune } 772ab8d36c9SPeter Brune 773ab8d36c9SPeter Brune /*@ 774ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the 775ab8d36c9SPeter Brune interpolation from l-1 to the lth level 776ab8d36c9SPeter Brune 777ab8d36c9SPeter Brune Input Parameters: 778f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 779ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 780ab8d36c9SPeter Brune 781f6dfbefdSBarry Smith Output Parameter: 782ab8d36c9SPeter Brune . mat - the interpolation operator 783ab8d36c9SPeter Brune 784ab8d36c9SPeter Brune Level: advanced 785ab8d36c9SPeter Brune 786f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASGetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetRScale()` 787ab8d36c9SPeter Brune @*/ 788d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) 789d71ae5a4SJacob Faibussowitsch { 79022d28d08SBarry Smith SNES_FAS *fas; 791ab8d36c9SPeter Brune SNES levelsnes; 79222d28d08SBarry Smith 793ab8d36c9SPeter Brune PetscFunctionBegin; 794f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 795f833ba53SLisandro Dalcin PetscValidPointer(mat, 3); 7969566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 797ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 798ab8d36c9SPeter Brune *mat = fas->interpolate; 799ab8d36c9SPeter Brune PetscFunctionReturn(0); 800ab8d36c9SPeter Brune } 801ab8d36c9SPeter Brune 802ab8d36c9SPeter Brune /*@ 803f6dfbefdSBarry Smith SNESFASSetRestriction - Sets the matrix to be used to restrict the defect 804ab8d36c9SPeter Brune from level l to l-1. 805ab8d36c9SPeter Brune 806ab8d36c9SPeter Brune Input Parameters: 807f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 808ab8d36c9SPeter Brune . mat - the restriction matrix 809ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 810ab8d36c9SPeter Brune 811ab8d36c9SPeter Brune Level: advanced 812ab8d36c9SPeter Brune 813ab8d36c9SPeter Brune Notes: 814ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation 815ab8d36c9SPeter Brune for the same level. 816ab8d36c9SPeter Brune 817ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 818ab8d36c9SPeter Brune out from the matrix size which one it is. 819ab8d36c9SPeter Brune 820ab8d36c9SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 821ab8d36c9SPeter Brune is used. 822ab8d36c9SPeter Brune 823f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetInjection()` 824ab8d36c9SPeter Brune @*/ 825d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) 826d71ae5a4SJacob Faibussowitsch { 82722d28d08SBarry Smith SNES_FAS *fas; 828ab8d36c9SPeter Brune SNES levelsnes; 82922d28d08SBarry Smith 830ab8d36c9SPeter Brune PetscFunctionBegin; 831f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 832f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 8339566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 834ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 8359566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 8369566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->restrct)); 837ab8d36c9SPeter Brune fas->restrct = mat; 838ab8d36c9SPeter Brune PetscFunctionReturn(0); 839ab8d36c9SPeter Brune } 840ab8d36c9SPeter Brune 841ab8d36c9SPeter Brune /*@ 842ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the 843ab8d36c9SPeter Brune restriction from l to the l-1th level 844ab8d36c9SPeter Brune 845ab8d36c9SPeter Brune Input Parameters: 846f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 847ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 848ab8d36c9SPeter Brune 849f6dfbefdSBarry Smith Output Parameter: 850ab8d36c9SPeter Brune . mat - the interpolation operator 851ab8d36c9SPeter Brune 852ab8d36c9SPeter Brune Level: advanced 853ab8d36c9SPeter Brune 854f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetRestriction()`, `SNESFASGetInjection()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()` 855ab8d36c9SPeter Brune @*/ 856d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) 857d71ae5a4SJacob Faibussowitsch { 85822d28d08SBarry Smith SNES_FAS *fas; 859ab8d36c9SPeter Brune SNES levelsnes; 86022d28d08SBarry Smith 861ab8d36c9SPeter Brune PetscFunctionBegin; 862f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 863f833ba53SLisandro Dalcin PetscValidPointer(mat, 3); 8649566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 865ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 866ab8d36c9SPeter Brune *mat = fas->restrct; 867ab8d36c9SPeter Brune PetscFunctionReturn(0); 868ab8d36c9SPeter Brune } 869ab8d36c9SPeter Brune 870ab8d36c9SPeter Brune /*@ 871ab8d36c9SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 872ab8d36c9SPeter Brune from level l to l-1. 873ab8d36c9SPeter Brune 874ab8d36c9SPeter Brune Input Parameters: 875f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 876ab8d36c9SPeter Brune . mat - the restriction matrix 877ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 878ab8d36c9SPeter Brune 879ab8d36c9SPeter Brune Level: advanced 880ab8d36c9SPeter Brune 881f6dfbefdSBarry Smith Note: 882ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to 883ab8d36c9SPeter Brune project the solution instead. 884ab8d36c9SPeter Brune 885f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetRestriction()` 886ab8d36c9SPeter Brune @*/ 887d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) 888d71ae5a4SJacob Faibussowitsch { 88922d28d08SBarry Smith SNES_FAS *fas; 890ab8d36c9SPeter Brune SNES levelsnes; 89122d28d08SBarry Smith 892ab8d36c9SPeter Brune PetscFunctionBegin; 893f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 894f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 8959566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 896ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 8979566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 8989566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->inject)); 8991aa26658SKarl Rupp 900ab8d36c9SPeter Brune fas->inject = mat; 901ab8d36c9SPeter Brune PetscFunctionReturn(0); 902ab8d36c9SPeter Brune } 903ab8d36c9SPeter Brune 904ab8d36c9SPeter Brune /*@ 905ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the 906ab8d36c9SPeter Brune injection from l-1 to the lth level 907ab8d36c9SPeter Brune 908ab8d36c9SPeter Brune Input Parameters: 909f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 910ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 911ab8d36c9SPeter Brune 912f6dfbefdSBarry Smith Output Parameter: 913ab8d36c9SPeter Brune . mat - the injection operator 914ab8d36c9SPeter Brune 915ab8d36c9SPeter Brune Level: advanced 916ab8d36c9SPeter Brune 917f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()` 918ab8d36c9SPeter Brune @*/ 919d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) 920d71ae5a4SJacob Faibussowitsch { 92122d28d08SBarry Smith SNES_FAS *fas; 922ab8d36c9SPeter Brune SNES levelsnes; 92322d28d08SBarry Smith 924ab8d36c9SPeter Brune PetscFunctionBegin; 925f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 926f833ba53SLisandro Dalcin PetscValidPointer(mat, 3); 9279566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 928ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 929ab8d36c9SPeter Brune *mat = fas->inject; 930ab8d36c9SPeter Brune PetscFunctionReturn(0); 931ab8d36c9SPeter Brune } 932ab8d36c9SPeter Brune 933ab8d36c9SPeter Brune /*@ 934ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 935ab8d36c9SPeter Brune operator from level l to l-1. 936ab8d36c9SPeter Brune 937ab8d36c9SPeter Brune Input Parameters: 938f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 939ab8d36c9SPeter Brune . rscale - the restriction scaling 940ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 941ab8d36c9SPeter Brune 942ab8d36c9SPeter Brune Level: advanced 943ab8d36c9SPeter Brune 944f6dfbefdSBarry Smith Note: 945ab8d36c9SPeter Brune This is only used in the case that the injection is not set. 946ab8d36c9SPeter Brune 947f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 948ab8d36c9SPeter Brune @*/ 949d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) 950d71ae5a4SJacob Faibussowitsch { 951ab8d36c9SPeter Brune SNES_FAS *fas; 952ab8d36c9SPeter Brune SNES levelsnes; 95322d28d08SBarry Smith 954ab8d36c9SPeter Brune PetscFunctionBegin; 955f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 956f833ba53SLisandro Dalcin if (rscale) PetscValidHeaderSpecific(rscale, VEC_CLASSID, 3); 9579566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 958ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 9599566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rscale)); 9609566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->rscale)); 961ab8d36c9SPeter Brune fas->rscale = rscale; 962ab8d36c9SPeter Brune PetscFunctionReturn(0); 963ab8d36c9SPeter Brune } 964ab8d36c9SPeter Brune 965ab8d36c9SPeter Brune /*@ 966ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level. 967ab8d36c9SPeter Brune 968ab8d36c9SPeter Brune Input Parameters: 969f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 970ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 971ab8d36c9SPeter Brune 972f6dfbefdSBarry Smith Output Parameter: 973ab8d36c9SPeter Brune smooth - the smoother 974ab8d36c9SPeter Brune 975ab8d36c9SPeter Brune Level: advanced 976ab8d36c9SPeter Brune 977f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 978ab8d36c9SPeter Brune @*/ 979d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) 980d71ae5a4SJacob Faibussowitsch { 981ab8d36c9SPeter Brune SNES_FAS *fas; 982ab8d36c9SPeter Brune SNES levelsnes; 98322d28d08SBarry Smith 984ab8d36c9SPeter Brune PetscFunctionBegin; 985f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 986f833ba53SLisandro Dalcin PetscValidPointer(smooth, 3); 9879566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 988ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 98948a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 990ab8d36c9SPeter Brune *smooth = fas->smoothd; 991ab8d36c9SPeter Brune PetscFunctionReturn(0); 992ab8d36c9SPeter Brune } 993ab8d36c9SPeter Brune 994ab8d36c9SPeter Brune /*@ 995ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level. 996ab8d36c9SPeter Brune 997ab8d36c9SPeter Brune Input Parameters: 998f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 999ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1000ab8d36c9SPeter Brune 1001f6dfbefdSBarry Smith Output Parameter: 1002ab8d36c9SPeter Brune smooth - the smoother 1003ab8d36c9SPeter Brune 1004ab8d36c9SPeter Brune Level: advanced 1005ab8d36c9SPeter Brune 1006f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1007ab8d36c9SPeter Brune @*/ 1008d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) 1009d71ae5a4SJacob Faibussowitsch { 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->smoothd; 1022ab8d36c9SPeter Brune PetscFunctionReturn(0); 1023ab8d36c9SPeter Brune } 1024ab8d36c9SPeter Brune 1025ab8d36c9SPeter Brune /*@ 1026ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level. 1027ab8d36c9SPeter Brune 1028ab8d36c9SPeter Brune Input Parameters: 1029f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1030ab8d36c9SPeter Brune - level - the level (0 is coarsest) 1031ab8d36c9SPeter Brune 1032f6dfbefdSBarry Smith Output Parameter: 1033ab8d36c9SPeter Brune smooth - the smoother 1034ab8d36c9SPeter Brune 1035ab8d36c9SPeter Brune Level: advanced 1036ab8d36c9SPeter Brune 1037f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1038ab8d36c9SPeter Brune @*/ 1039d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) 1040d71ae5a4SJacob Faibussowitsch { 1041ab8d36c9SPeter Brune SNES_FAS *fas; 1042ab8d36c9SPeter Brune SNES levelsnes; 104322d28d08SBarry Smith 1044ab8d36c9SPeter Brune PetscFunctionBegin; 1045f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1046f833ba53SLisandro Dalcin PetscValidPointer(smooth, 3); 10479566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 1048ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1049ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 105048a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 105148a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu)); 1052ab8d36c9SPeter Brune *smooth = fas->smoothu; 1053ab8d36c9SPeter Brune PetscFunctionReturn(0); 1054ab8d36c9SPeter Brune } 1055d6ad1212SPeter Brune 1056d6ad1212SPeter Brune /*@ 1057d6ad1212SPeter Brune SNESFASGetCoarseSolve - Gets the coarsest solver. 1058d6ad1212SPeter Brune 1059f6dfbefdSBarry Smith Input Parameter: 1060f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 1061d6ad1212SPeter Brune 1062f6dfbefdSBarry Smith Output Parameter: 1063a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver 1064d6ad1212SPeter Brune 1065d6ad1212SPeter Brune Level: advanced 1066d6ad1212SPeter Brune 1067f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1068d6ad1212SPeter Brune @*/ 1069d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse) 1070d71ae5a4SJacob Faibussowitsch { 1071d6ad1212SPeter Brune SNES_FAS *fas; 1072d6ad1212SPeter Brune SNES levelsnes; 107322d28d08SBarry Smith 1074d6ad1212SPeter Brune PetscFunctionBegin; 1075f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1076064a246eSJacob Faibussowitsch PetscValidPointer(coarse, 2); 10779566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, 0, &levelsnes)); 1078d6ad1212SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1079d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 108048a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 1081a3a80b83SMatthew G. Knepley *coarse = fas->smoothd; 1082d6ad1212SPeter Brune PetscFunctionReturn(0); 1083d6ad1212SPeter Brune } 1084928e959bSPeter Brune 1085928e959bSPeter Brune /*@ 1086f6dfbefdSBarry Smith SNESFASFullSetDownSweep - Smooth during the initial downsweep for `SNESFAS` 1087928e959bSPeter Brune 1088*c3339decSBarry Smith Logically Collective 1089928e959bSPeter Brune 1090928e959bSPeter Brune Input Parameters: 1091f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1092928e959bSPeter Brune - swp - whether to downsweep or not 1093928e959bSPeter Brune 1094928e959bSPeter Brune Options Database Key: 1095928e959bSPeter Brune . -snes_fas_full_downsweep - Sets number of pre-smoothing steps 1096928e959bSPeter Brune 1097928e959bSPeter Brune Level: advanced 1098928e959bSPeter Brune 1099f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 1100928e959bSPeter Brune @*/ 1101d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullSetDownSweep(SNES snes, PetscBool swp) 1102d71ae5a4SJacob Faibussowitsch { 1103f833ba53SLisandro Dalcin SNES_FAS *fas; 1104928e959bSPeter Brune 1105928e959bSPeter Brune PetscFunctionBegin; 1106f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1107f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 1108928e959bSPeter Brune fas->full_downsweep = swp; 11091baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetDownSweep(fas->next, swp)); 1110928e959bSPeter Brune PetscFunctionReturn(0); 1111928e959bSPeter Brune } 11126dfbafefSToby Isaac 11136dfbafefSToby Isaac /*@ 11146dfbafefSToby Isaac SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 11156dfbafefSToby Isaac 1116*c3339decSBarry Smith Logically Collective 11176dfbafefSToby Isaac 11186dfbafefSToby Isaac Input Parameters: 1119f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 11206dfbafefSToby Isaac - total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 11216dfbafefSToby Isaac 11226dfbafefSToby Isaac Options Database Key: 112367b8a455SSatish Balay . -snes_fas_full_total - Use total restriction and interpolation on the initial down and up sweeps for the full FAS cycle 11246dfbafefSToby Isaac 11256dfbafefSToby Isaac Level: advanced 11266dfbafefSToby Isaac 1127f6dfbefdSBarry Smith Note: 1128f6dfbefdSBarry Smith This option is only significant if the interpolation of a coarse correction (`MatInterpolate()`) is significantly different from total 1129f6dfbefdSBarry Smith solution interpolation (`DMInterpolateSolution()`). 11306dfbafefSToby Isaac 1131f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()` 11326dfbafefSToby Isaac @*/ 1133d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullSetTotal(SNES snes, PetscBool total) 1134d71ae5a4SJacob Faibussowitsch { 11356dfbafefSToby Isaac SNES_FAS *fas; 11366dfbafefSToby Isaac 11376dfbafefSToby Isaac PetscFunctionBegin; 11386dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 11396dfbafefSToby Isaac fas = (SNES_FAS *)snes->data; 11406dfbafefSToby Isaac fas->full_total = total; 11411baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetTotal(fas->next, total)); 11426dfbafefSToby Isaac PetscFunctionReturn(0); 11436dfbafefSToby Isaac } 11446dfbafefSToby Isaac 11456dfbafefSToby Isaac /*@ 11466dfbafefSToby Isaac SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 11476dfbafefSToby Isaac 1148*c3339decSBarry Smith Logically Collective 11496dfbafefSToby Isaac 1150f6dfbefdSBarry Smith Input Parameter: 1151f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 11526dfbafefSToby Isaac 11536dfbafefSToby Isaac Output: 11546dfbafefSToby Isaac . total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 11556dfbafefSToby Isaac 11566dfbafefSToby Isaac Level: advanced 11576dfbafefSToby Isaac 1158f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`, `SNESFullSetTotal()` 11596dfbafefSToby Isaac @*/ 1160d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullGetTotal(SNES snes, PetscBool *total) 1161d71ae5a4SJacob Faibussowitsch { 11626dfbafefSToby Isaac SNES_FAS *fas; 11636dfbafefSToby Isaac 11646dfbafefSToby Isaac PetscFunctionBegin; 11656dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 11666dfbafefSToby Isaac fas = (SNES_FAS *)snes->data; 11676dfbafefSToby Isaac *total = fas->full_total; 11686dfbafefSToby Isaac PetscFunctionReturn(0); 11696dfbafefSToby Isaac } 1170