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)); 263ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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); 504f572ea9SToby Isaac PetscAssertPointer(fastype, 2); 51f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 52ab8d36c9SPeter Brune *fastype = fas->fastype; 533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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: 692fe279fdSBarry Smith If the number of levels is one then the multigrid uses the `-fas_levels` prefix 702fe279fdSBarry Smith 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) { 883ba16761SJacob Faibussowitsch if (!comms) PetscFunctionReturn(PETSC_SUCCESS); 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 } 1233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 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 132e4094ef1SJacob Faibussowitsch 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); 1454f572ea9SToby Isaac PetscAssertPointer(levels, 2); 146f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 147ab8d36c9SPeter Brune *levels = fas->levels; 1483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 149ab8d36c9SPeter Brune } 150ab8d36c9SPeter Brune 151ab8d36c9SPeter Brune /*@ 152*ceaaa498SBarry Smith SNESFASGetCycleSNES - Gets the `SNES` corresponding to a particular level of the `SNESFAS` hierarchy 153ab8d36c9SPeter Brune 154ab8d36c9SPeter Brune Input Parameters: 155f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context 156f6dfbefdSBarry Smith - level - the level to get 157f6dfbefdSBarry Smith 158f6dfbefdSBarry Smith Output Parameter: 159f6dfbefdSBarry Smith . lsnes - the `SNES` for the requested level 160ab8d36c9SPeter Brune 161ab8d36c9SPeter Brune Level: advanced 162ab8d36c9SPeter Brune 163f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()`, `SNESFASGetLevels()` 164ab8d36c9SPeter Brune @*/ 165d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetCycleSNES(SNES snes, PetscInt level, SNES *lsnes) 166d71ae5a4SJacob Faibussowitsch { 167f833ba53SLisandro Dalcin SNES_FAS *fas; 168ab8d36c9SPeter Brune PetscInt i; 169ab8d36c9SPeter Brune 170ab8d36c9SPeter Brune PetscFunctionBegin; 171f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1724f572ea9SToby Isaac PetscAssertPointer(lsnes, 3); 173f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 17463a3b9bcSJacob 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); 1750b121fc5SBarry Smith PetscCheck(fas->level == fas->levels - 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetCycleSNES may only be called on the finest-level SNES"); 176ab8d36c9SPeter Brune 177ab8d36c9SPeter Brune *lsnes = snes; 178ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) { 179ab8d36c9SPeter Brune *lsnes = fas->next; 180ab8d36c9SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 181ab8d36c9SPeter Brune } 18208401ef6SPierre Jolivet PetscCheck(fas->level == level, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "SNESFAS level hierarchy corrupt"); 1833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 184ab8d36c9SPeter Brune } 185ab8d36c9SPeter Brune 186ab8d36c9SPeter Brune /*@ 187ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 188ab8d36c9SPeter Brune use on all levels. 189ab8d36c9SPeter Brune 190c3339decSBarry Smith Logically Collective 191ab8d36c9SPeter Brune 192ab8d36c9SPeter Brune Input Parameters: 193f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context 194ab8d36c9SPeter Brune - n - the number of smoothing steps 195ab8d36c9SPeter Brune 196ab8d36c9SPeter Brune Options Database Key: 197ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 198ab8d36c9SPeter Brune 199ab8d36c9SPeter Brune Level: advanced 200ab8d36c9SPeter Brune 201f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothDown()` 202ab8d36c9SPeter Brune @*/ 203d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) 204d71ae5a4SJacob Faibussowitsch { 205f833ba53SLisandro Dalcin SNES_FAS *fas; 20622d28d08SBarry Smith 207ab8d36c9SPeter Brune PetscFunctionBegin; 208f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 209f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 210ab8d36c9SPeter Brune fas->max_up_it = n; 21148a46eb9SPierre Jolivet if (!fas->smoothu && fas->level != 0) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu)); 2121baa6e33SBarry Smith if (fas->smoothu) PetscCall(SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs)); 2131baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothUp(fas->next, n)); 2143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 215ab8d36c9SPeter Brune } 216ab8d36c9SPeter Brune 217ab8d36c9SPeter Brune /*@ 218ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 219ab8d36c9SPeter Brune use on all levels. 220ab8d36c9SPeter Brune 221c3339decSBarry Smith Logically Collective 222ab8d36c9SPeter Brune 223ab8d36c9SPeter Brune Input Parameters: 224f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 225ab8d36c9SPeter Brune - n - the number of smoothing steps 226ab8d36c9SPeter Brune 227ab8d36c9SPeter Brune Options Database Key: 228ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 229ab8d36c9SPeter Brune 230ab8d36c9SPeter Brune Level: advanced 231ab8d36c9SPeter Brune 232f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 233ab8d36c9SPeter Brune @*/ 234d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) 235d71ae5a4SJacob Faibussowitsch { 236f833ba53SLisandro Dalcin SNES_FAS *fas; 23722d28d08SBarry Smith 238ab8d36c9SPeter Brune PetscFunctionBegin; 239f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 240f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 24148a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd)); 2429566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs)); 2431aa26658SKarl Rupp 244ab8d36c9SPeter Brune fas->max_down_it = n; 2451baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothDown(fas->next, n)); 2463ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 247ab8d36c9SPeter Brune } 248ab8d36c9SPeter Brune 24987f44e3fSPeter Brune /*@ 250f6dfbefdSBarry Smith SNESFASSetContinuation - Sets the `SNESFAS` cycle to default to exact Newton solves on the upsweep 25187f44e3fSPeter Brune 252c3339decSBarry Smith Logically Collective 25387f44e3fSPeter Brune 25487f44e3fSPeter Brune Input Parameters: 255f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 256e4094ef1SJacob Faibussowitsch - continuation - the number of smoothing steps 25787f44e3fSPeter Brune 25887f44e3fSPeter Brune Options Database Key: 25987f44e3fSPeter Brune . -snes_fas_continuation - sets continuation to true 26087f44e3fSPeter Brune 26187f44e3fSPeter Brune Level: advanced 26287f44e3fSPeter Brune 263f6dfbefdSBarry Smith Note: 26495452b02SPatrick Sanan This sets the prefix on the upsweep smoothers to -fas_continuation 26587f44e3fSPeter Brune 266f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 26787f44e3fSPeter Brune @*/ 268d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetContinuation(SNES snes, PetscBool continuation) 269d71ae5a4SJacob Faibussowitsch { 27087f44e3fSPeter Brune const char *optionsprefix; 27187f44e3fSPeter Brune char tprefix[128]; 272f833ba53SLisandro Dalcin SNES_FAS *fas; 27387f44e3fSPeter Brune 27487f44e3fSPeter Brune PetscFunctionBegin; 275f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 276f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 2779566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 27848a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu)); 2799566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_levels_continuation_", sizeof(tprefix))); 2809566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(fas->smoothu, optionsprefix)); 2819566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->smoothu, tprefix)); 2829566063dSJacob Faibussowitsch PetscCall(SNESSetType(fas->smoothu, SNESNEWTONLS)); 2839566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothu, fas->fine->abstol, fas->fine->rtol, fas->fine->stol, 50, 100)); 28487f44e3fSPeter Brune fas->continuation = continuation; 2851baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetContinuation(fas->next, continuation)); 2863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28787f44e3fSPeter Brune } 28887f44e3fSPeter Brune 289ab8d36c9SPeter Brune /*@ 290f6dfbefdSBarry Smith SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use `SNESFASSetCyclesOnLevel()` for more 291ab8d36c9SPeter Brune complicated cycling. 292ab8d36c9SPeter Brune 293c3339decSBarry Smith Logically Collective 294ab8d36c9SPeter Brune 295ab8d36c9SPeter Brune Input Parameters: 296f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 297ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 298ab8d36c9SPeter Brune 299ab8d36c9SPeter Brune Options Database Key: 30067b8a455SSatish Balay . -snes_fas_cycles <1,2> - 1 for V-cycle, 2 for W-cycle 301ab8d36c9SPeter Brune 302ab8d36c9SPeter Brune Level: advanced 303ab8d36c9SPeter Brune 304f6dfbefdSBarry Smith .seealso: `SNES`, `SNESFAS`, `SNESFASSetCyclesOnLevel()` 305ab8d36c9SPeter Brune @*/ 306d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) 307d71ae5a4SJacob Faibussowitsch { 308f833ba53SLisandro Dalcin SNES_FAS *fas; 309ab8d36c9SPeter Brune PetscBool isFine; 31022d28d08SBarry Smith 311ab8d36c9SPeter Brune PetscFunctionBegin; 312f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3139566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 314f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 315ab8d36c9SPeter Brune fas->n_cycles = cycles; 31648a46eb9SPierre Jolivet if (!isFine) PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs)); 3171baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetCycles(fas->next, cycles)); 3183ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 319ab8d36c9SPeter Brune } 320ab8d36c9SPeter Brune 321c8c899caSPeter Brune /*@ 322c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring 323c8c899caSPeter Brune 324c3339decSBarry Smith Logically Collective 325c8c899caSPeter Brune 326c8c899caSPeter Brune Input Parameters: 327f6dfbefdSBarry Smith + snes - the `SNESFAS` context 3282fe279fdSBarry Smith . vf - viewer and format structure (may be `NULL` if flg is `PETSC_FALSE`) 329c8c899caSPeter Brune - flg - monitor or not 330c8c899caSPeter Brune 331c8c899caSPeter Brune Level: advanced 332c8c899caSPeter Brune 333f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESSetMonitor()`, `SNESFASSetCyclesOnLevel()` 334c8c899caSPeter Brune @*/ 335d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg) 336d71ae5a4SJacob Faibussowitsch { 337f833ba53SLisandro Dalcin SNES_FAS *fas; 338c8c899caSPeter Brune PetscBool isFine; 339f833ba53SLisandro Dalcin PetscInt i, levels; 340c8c899caSPeter Brune SNES levelsnes; 34122d28d08SBarry Smith 342c8c899caSPeter Brune PetscFunctionBegin; 343f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3449566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 345f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 346f833ba53SLisandro Dalcin levels = fas->levels; 347c8c899caSPeter Brune if (isFine) { 348c8c899caSPeter Brune for (i = 0; i < levels; i++) { 3499566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 350c8c899caSPeter Brune fas = (SNES_FAS *)levelsnes->data; 351c8c899caSPeter Brune if (flg) { 352c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */ 3539566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes)); 354d142ab34SLawrence Mitchell /* Only register destroy on finest level */ 3559566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(levelsnes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))SNESMonitorDefault, vf, (!i ? (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy : NULL))); 3561aa26658SKarl Rupp } else if (i != fas->levels - 1) { 357c8c899caSPeter Brune /* unset the monitors on the coarse levels */ 3589566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes)); 359c8c899caSPeter Brune } 360c8c899caSPeter Brune } 361c8c899caSPeter Brune } 3623ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 363c8c899caSPeter Brune } 364c8c899caSPeter Brune 3650dd27c6cSPeter Brune /*@ 366f6dfbefdSBarry Smith SNESFASSetLog - Sets or unsets time logging for various `SNESFAS` stages on all levels 3670dd27c6cSPeter Brune 368c3339decSBarry Smith Logically Collective 3690dd27c6cSPeter Brune 3700dd27c6cSPeter Brune Input Parameters: 371f6dfbefdSBarry Smith + snes - the `SNESFAS` context 3720dd27c6cSPeter Brune - flg - monitor or not 3730dd27c6cSPeter Brune 3740dd27c6cSPeter Brune Level: advanced 3750dd27c6cSPeter Brune 376f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetMonitor()` 3770dd27c6cSPeter Brune @*/ 378d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) 379d71ae5a4SJacob Faibussowitsch { 380f833ba53SLisandro Dalcin SNES_FAS *fas; 3810dd27c6cSPeter Brune PetscBool isFine; 382f833ba53SLisandro Dalcin PetscInt i, levels; 3830dd27c6cSPeter Brune SNES levelsnes; 3840dd27c6cSPeter Brune char eventname[128]; 3850dd27c6cSPeter Brune 3860dd27c6cSPeter Brune PetscFunctionBegin; 387f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3889566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 389f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 390f833ba53SLisandro Dalcin levels = fas->levels; 3910dd27c6cSPeter Brune if (isFine) { 3920dd27c6cSPeter Brune for (i = 0; i < levels; i++) { 3939566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 3940dd27c6cSPeter Brune fas = (SNES_FAS *)levelsnes->data; 3950dd27c6cSPeter Brune if (flg) { 3969566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSetup %d", (int)i)); 3979566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsetup)); 3989566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSmooth %d", (int)i)); 3999566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsolve)); 4009566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASResid %d", (int)i)); 4019566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventresidual)); 4029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASInterp %d", (int)i)); 4039566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventinterprestrict)); 4040dd27c6cSPeter Brune } else { 4050298fd71SBarry Smith fas->eventsmoothsetup = 0; 4060298fd71SBarry Smith fas->eventsmoothsolve = 0; 4070298fd71SBarry Smith fas->eventresidual = 0; 4080298fd71SBarry Smith fas->eventinterprestrict = 0; 4090dd27c6cSPeter Brune } 4100dd27c6cSPeter Brune } 4110dd27c6cSPeter Brune } 4123ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4130dd27c6cSPeter Brune } 4140dd27c6cSPeter Brune 415ab8d36c9SPeter Brune /* 416ab8d36c9SPeter Brune Creates the default smoother type. 417ab8d36c9SPeter Brune 41804d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level. 419ab8d36c9SPeter Brune 420ab8d36c9SPeter Brune */ 421d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) 422d71ae5a4SJacob Faibussowitsch { 423ab8d36c9SPeter Brune SNES_FAS *fas; 424ab8d36c9SPeter Brune const char *optionsprefix; 425ab8d36c9SPeter Brune char tprefix[128]; 426ab8d36c9SPeter Brune SNES nsmooth; 42722d28d08SBarry Smith 428ab8d36c9SPeter Brune PetscFunctionBegin; 429f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 4304f572ea9SToby Isaac PetscAssertPointer(smooth, 2); 431ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 4329566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 433ab8d36c9SPeter Brune /* create the default smoother */ 4349566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth)); 435ab8d36c9SPeter Brune if (fas->level == 0) { 4369566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_coarse_", sizeof(tprefix))); 4379566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix)); 4389566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix)); 4399566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNEWTONLS)); 4409566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs)); 441ab8d36c9SPeter Brune } else { 4429566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%d_", (int)fas->level)); 4439566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix)); 4449566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix)); 4459566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNRICHARDSON)); 4469566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs)); 447ab8d36c9SPeter Brune } 4489566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1)); 4499566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth)); 4509566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)nsmooth, PetscMGLevelId, fas->level)); 451ab8d36c9SPeter Brune *smooth = nsmooth; 4523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 453ab8d36c9SPeter Brune } 454ab8d36c9SPeter Brune 455ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */ 456ab8d36c9SPeter Brune 457ab8d36c9SPeter Brune /*@ 458*ceaaa498SBarry Smith SNESFASCycleSetCycles - Sets the number of cycles for all levels in a `SNESFAS` 459ab8d36c9SPeter Brune 460c3339decSBarry Smith Logically Collective 461ab8d36c9SPeter Brune 462ab8d36c9SPeter Brune Input Parameters: 463f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 464ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 465ab8d36c9SPeter Brune 466ab8d36c9SPeter Brune Level: advanced 467ab8d36c9SPeter Brune 468f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetCycles()` 469ab8d36c9SPeter Brune @*/ 470d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) 471d71ae5a4SJacob Faibussowitsch { 472f833ba53SLisandro Dalcin SNES_FAS *fas; 47322d28d08SBarry Smith 474ab8d36c9SPeter Brune PetscFunctionBegin; 475f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 476f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 477ab8d36c9SPeter Brune fas->n_cycles = cycles; 4789566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs)); 4793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 480ab8d36c9SPeter Brune } 481ab8d36c9SPeter Brune 482ab8d36c9SPeter Brune /*@ 483ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 484ab8d36c9SPeter Brune 485c3339decSBarry Smith Logically Collective 486ab8d36c9SPeter Brune 487f6dfbefdSBarry Smith Input Parameter: 488f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 489ab8d36c9SPeter Brune 490f6dfbefdSBarry Smith Output Parameter: 491ab8d36c9SPeter Brune . smooth - the smoother 492ab8d36c9SPeter Brune 493ab8d36c9SPeter Brune Level: advanced 494ab8d36c9SPeter Brune 495*ceaaa498SBarry Smith Note: 496*ceaaa498SBarry Smith The `snes` should be obtained with `SNESFASGetCycleSNES()` 497*ceaaa498SBarry Smith 498*ceaaa498SBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()`, `SNESFASGetCycleSNES()` 499ab8d36c9SPeter Brune @*/ 500d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 501d71ae5a4SJacob Faibussowitsch { 502ab8d36c9SPeter Brune SNES_FAS *fas; 5035fd66863SKarl Rupp 504ab8d36c9SPeter Brune PetscFunctionBegin; 505f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 5064f572ea9SToby Isaac PetscAssertPointer(smooth, 2); 507ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 508ab8d36c9SPeter Brune *smooth = fas->smoothd; 5093ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 510ab8d36c9SPeter Brune } 511ab8d36c9SPeter Brune /*@ 512ab8d36c9SPeter Brune SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 513ab8d36c9SPeter Brune 514c3339decSBarry Smith Logically Collective 515ab8d36c9SPeter Brune 516f6dfbefdSBarry Smith Input Parameter: 517f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 518ab8d36c9SPeter Brune 519f6dfbefdSBarry Smith Output Parameter: 520ab8d36c9SPeter Brune . smoothu - the smoother 521ab8d36c9SPeter Brune 522*ceaaa498SBarry Smith Level: advanced 523*ceaaa498SBarry Smith 524f6dfbefdSBarry Smith Note: 525ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent 526ab8d36c9SPeter Brune default behavior in the process of the solve. 527ab8d36c9SPeter Brune 528*ceaaa498SBarry Smith The `snes` should be obtained with `SNESFASGetCycleSNES()` 529ab8d36c9SPeter Brune 530*ceaaa498SBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherDown()`, `SNESFASGetCycleSNES()` 531ab8d36c9SPeter Brune @*/ 532d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 533d71ae5a4SJacob Faibussowitsch { 534ab8d36c9SPeter Brune SNES_FAS *fas; 5355fd66863SKarl Rupp 536ab8d36c9SPeter Brune PetscFunctionBegin; 537f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 5384f572ea9SToby Isaac PetscAssertPointer(smoothu, 2); 539ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 5401aa26658SKarl Rupp if (!fas->smoothu) *smoothu = fas->smoothd; 5411aa26658SKarl Rupp else *smoothu = fas->smoothu; 5423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 543ab8d36c9SPeter Brune } 544ab8d36c9SPeter Brune 545ab8d36c9SPeter Brune /*@ 546ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 547ab8d36c9SPeter Brune 548c3339decSBarry Smith Logically Collective 549ab8d36c9SPeter Brune 550f6dfbefdSBarry Smith Input Parameter: 551f6dfbefdSBarry Smith . snes - `SNESFAS`, the nonlinear multigrid context 552ab8d36c9SPeter Brune 553f6dfbefdSBarry Smith Output Parameter: 554ab8d36c9SPeter Brune . smoothd - the smoother 555ab8d36c9SPeter Brune 556ab8d36c9SPeter Brune Level: advanced 557ab8d36c9SPeter Brune 558*ceaaa498SBarry Smith Note: 559*ceaaa498SBarry Smith The `snes` should be obtained with `SNESFASGetCycleSNES()` 560*ceaaa498SBarry Smith 561*ceaaa498SBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`, `SNESFASGetCycleSNES()` 562ab8d36c9SPeter Brune @*/ 563d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 564d71ae5a4SJacob Faibussowitsch { 565ab8d36c9SPeter Brune SNES_FAS *fas; 5665fd66863SKarl Rupp 567ab8d36c9SPeter Brune PetscFunctionBegin; 568f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 5694f572ea9SToby Isaac PetscAssertPointer(smoothd, 2); 570ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 571ab8d36c9SPeter Brune *smoothd = fas->smoothd; 5723ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 573ab8d36c9SPeter Brune } 574ab8d36c9SPeter Brune 575ab8d36c9SPeter Brune /*@ 576ab8d36c9SPeter Brune SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 577ab8d36c9SPeter Brune 578c3339decSBarry Smith Logically Collective 579ab8d36c9SPeter Brune 580f6dfbefdSBarry Smith Input Parameter: 581f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 582ab8d36c9SPeter Brune 583f6dfbefdSBarry Smith Output Parameter: 584f6dfbefdSBarry Smith . correction - the coarse correction solve on this level 585ab8d36c9SPeter Brune 586f6dfbefdSBarry Smith Note: 5870298fd71SBarry Smith Returns NULL on the coarsest level. 588ab8d36c9SPeter Brune 589ab8d36c9SPeter Brune Level: advanced 590ab8d36c9SPeter Brune 591f6dfbefdSBarry Smith .seealso: `SNESFAS` `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 592ab8d36c9SPeter Brune @*/ 593d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 594d71ae5a4SJacob Faibussowitsch { 595ab8d36c9SPeter Brune SNES_FAS *fas; 5965fd66863SKarl Rupp 597ab8d36c9SPeter Brune PetscFunctionBegin; 598f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 5994f572ea9SToby Isaac PetscAssertPointer(correction, 2); 600ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 601ab8d36c9SPeter Brune *correction = fas->next; 6023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 603ab8d36c9SPeter Brune } 604ab8d36c9SPeter Brune 605ab8d36c9SPeter Brune /*@ 606ab8d36c9SPeter Brune SNESFASCycleGetInterpolation - Gets the interpolation on this level 607ab8d36c9SPeter Brune 608c3339decSBarry Smith Logically Collective 609ab8d36c9SPeter Brune 610f6dfbefdSBarry Smith Input Parameter: 611f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 612ab8d36c9SPeter Brune 613f6dfbefdSBarry Smith Output Parameter: 614ab8d36c9SPeter Brune . mat - the interpolation operator on this level 615ab8d36c9SPeter Brune 616f6dfbefdSBarry Smith Level: advanced 617ab8d36c9SPeter Brune 618f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 619ab8d36c9SPeter Brune @*/ 620d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 621d71ae5a4SJacob Faibussowitsch { 622ab8d36c9SPeter Brune SNES_FAS *fas; 6235fd66863SKarl Rupp 624ab8d36c9SPeter Brune PetscFunctionBegin; 625f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 6264f572ea9SToby Isaac PetscAssertPointer(mat, 2); 627ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 628ab8d36c9SPeter Brune *mat = fas->interpolate; 6293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 630ab8d36c9SPeter Brune } 631ab8d36c9SPeter Brune 632ab8d36c9SPeter Brune /*@ 633ab8d36c9SPeter Brune SNESFASCycleGetRestriction - Gets the restriction on this level 634ab8d36c9SPeter Brune 635c3339decSBarry Smith Logically Collective 636ab8d36c9SPeter Brune 637f6dfbefdSBarry Smith Input Parameter: 638f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 639ab8d36c9SPeter Brune 640f6dfbefdSBarry Smith Output Parameter: 641ab8d36c9SPeter Brune . mat - the restriction operator on this level 642ab8d36c9SPeter Brune 643f6dfbefdSBarry Smith Level: advanced 644ab8d36c9SPeter Brune 645f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetRestriction()`, `SNESFASCycleGetInterpolation()` 646ab8d36c9SPeter Brune @*/ 647d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 648d71ae5a4SJacob Faibussowitsch { 649ab8d36c9SPeter Brune SNES_FAS *fas; 6505fd66863SKarl Rupp 651ab8d36c9SPeter Brune PetscFunctionBegin; 652f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 6534f572ea9SToby Isaac PetscAssertPointer(mat, 2); 654ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 655ab8d36c9SPeter Brune *mat = fas->restrct; 6563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 657ab8d36c9SPeter Brune } 658ab8d36c9SPeter Brune 659ab8d36c9SPeter Brune /*@ 660ab8d36c9SPeter Brune SNESFASCycleGetInjection - Gets the injection on this level 661ab8d36c9SPeter Brune 662c3339decSBarry Smith Logically Collective 663ab8d36c9SPeter Brune 664f6dfbefdSBarry Smith Input Parameter: 665f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 666ab8d36c9SPeter Brune 667f6dfbefdSBarry Smith Output Parameter: 668ab8d36c9SPeter Brune . mat - the restriction operator on this level 669ab8d36c9SPeter Brune 670f6dfbefdSBarry Smith Level: advanced 671ab8d36c9SPeter Brune 672f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetInjection()`, `SNESFASCycleGetRestriction()` 673ab8d36c9SPeter Brune @*/ 674d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 675d71ae5a4SJacob Faibussowitsch { 676ab8d36c9SPeter Brune SNES_FAS *fas; 6775fd66863SKarl Rupp 678ab8d36c9SPeter Brune PetscFunctionBegin; 679f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 6804f572ea9SToby Isaac PetscAssertPointer(mat, 2); 681ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 682ab8d36c9SPeter Brune *mat = fas->inject; 6833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 684ab8d36c9SPeter Brune } 685ab8d36c9SPeter Brune 686ab8d36c9SPeter Brune /*@ 687ab8d36c9SPeter Brune SNESFASCycleGetRScale - Gets the injection on this level 688ab8d36c9SPeter Brune 689c3339decSBarry Smith Logically Collective 690ab8d36c9SPeter Brune 691f6dfbefdSBarry Smith Input Parameter: 692f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 693ab8d36c9SPeter Brune 694f6dfbefdSBarry Smith Output Parameter: 695e4094ef1SJacob Faibussowitsch . vec - the restriction operator on this level 696ab8d36c9SPeter Brune 697f6dfbefdSBarry Smith Level: advanced 698ab8d36c9SPeter Brune 699f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetRestriction()`, `SNESFASGetRScale()` 700ab8d36c9SPeter Brune @*/ 701d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 702d71ae5a4SJacob Faibussowitsch { 703ab8d36c9SPeter Brune SNES_FAS *fas; 7045fd66863SKarl Rupp 705ab8d36c9SPeter Brune PetscFunctionBegin; 706f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 7074f572ea9SToby Isaac PetscAssertPointer(vec, 2); 708ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 709ab8d36c9SPeter Brune *vec = fas->rscale; 7103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 711ab8d36c9SPeter Brune } 712ab8d36c9SPeter Brune 713ab8d36c9SPeter Brune /*@ 714ab8d36c9SPeter Brune SNESFASCycleIsFine - Determines if a given cycle is the fine level. 715ab8d36c9SPeter Brune 716c3339decSBarry Smith Logically Collective 717ab8d36c9SPeter Brune 718f6dfbefdSBarry Smith Input Parameter: 719f6dfbefdSBarry Smith . snes - the `SNESFAS` `SNES` context 720ab8d36c9SPeter Brune 721f6dfbefdSBarry Smith Output Parameter: 722ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not 723ab8d36c9SPeter Brune 724ab8d36c9SPeter Brune Level: advanced 725ab8d36c9SPeter Brune 726f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()` 727ab8d36c9SPeter Brune @*/ 728d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 729d71ae5a4SJacob Faibussowitsch { 730ab8d36c9SPeter Brune SNES_FAS *fas; 7315fd66863SKarl Rupp 732ab8d36c9SPeter Brune PetscFunctionBegin; 733f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 7344f572ea9SToby Isaac PetscAssertPointer(flg, 2); 735ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 7361aa26658SKarl Rupp if (fas->level == fas->levels - 1) *flg = PETSC_TRUE; 7371aa26658SKarl Rupp else *flg = PETSC_FALSE; 7383ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 739ab8d36c9SPeter Brune } 740ab8d36c9SPeter Brune 741f6dfbefdSBarry Smith /* functions called on the finest level that return level-specific information */ 742ab8d36c9SPeter Brune 743ab8d36c9SPeter Brune /*@ 744f6dfbefdSBarry Smith SNESFASSetInterpolation - Sets the `Mat` to be used to apply the 745ab8d36c9SPeter Brune interpolation from l-1 to the lth level 746ab8d36c9SPeter Brune 747ab8d36c9SPeter Brune Input Parameters: 748f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 749ab8d36c9SPeter Brune . mat - the interpolation operator 750ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 751ab8d36c9SPeter Brune 752ab8d36c9SPeter Brune Level: advanced 753ab8d36c9SPeter Brune 754ab8d36c9SPeter Brune Notes: 755ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction 756ab8d36c9SPeter Brune for the same level. 757ab8d36c9SPeter Brune 758ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 759ab8d36c9SPeter Brune out from the matrix size which one it is. 760ab8d36c9SPeter Brune 761f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`, `SNESFASSetRScale()` 762ab8d36c9SPeter Brune @*/ 763d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) 764d71ae5a4SJacob Faibussowitsch { 76522d28d08SBarry Smith SNES_FAS *fas; 766ab8d36c9SPeter Brune SNES levelsnes; 76722d28d08SBarry Smith 768ab8d36c9SPeter Brune PetscFunctionBegin; 769f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 770f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 7719566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 772ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 7739566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 7749566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->interpolate)); 775ab8d36c9SPeter Brune fas->interpolate = mat; 7763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 777ab8d36c9SPeter Brune } 778ab8d36c9SPeter Brune 779ab8d36c9SPeter Brune /*@ 780ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the 781ab8d36c9SPeter Brune interpolation from l-1 to the lth level 782ab8d36c9SPeter Brune 783ab8d36c9SPeter Brune Input Parameters: 784f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 785ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 786ab8d36c9SPeter Brune 787f6dfbefdSBarry Smith Output Parameter: 788ab8d36c9SPeter Brune . mat - the interpolation operator 789ab8d36c9SPeter Brune 790ab8d36c9SPeter Brune Level: advanced 791ab8d36c9SPeter Brune 792f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASGetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetRScale()` 793ab8d36c9SPeter Brune @*/ 794d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) 795d71ae5a4SJacob Faibussowitsch { 79622d28d08SBarry Smith SNES_FAS *fas; 797ab8d36c9SPeter Brune SNES levelsnes; 79822d28d08SBarry Smith 799ab8d36c9SPeter Brune PetscFunctionBegin; 800f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 8014f572ea9SToby Isaac PetscAssertPointer(mat, 3); 8029566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 803ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 804ab8d36c9SPeter Brune *mat = fas->interpolate; 8053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 806ab8d36c9SPeter Brune } 807ab8d36c9SPeter Brune 808ab8d36c9SPeter Brune /*@ 809f6dfbefdSBarry Smith SNESFASSetRestriction - Sets the matrix to be used to restrict the defect 810ab8d36c9SPeter Brune from level l to l-1. 811ab8d36c9SPeter Brune 812ab8d36c9SPeter Brune Input Parameters: 813f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 814ab8d36c9SPeter Brune . mat - the restriction matrix 815ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 816ab8d36c9SPeter Brune 817ab8d36c9SPeter Brune Level: advanced 818ab8d36c9SPeter Brune 819ab8d36c9SPeter Brune Notes: 820ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation 821ab8d36c9SPeter Brune for the same level. 822ab8d36c9SPeter Brune 823ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 824ab8d36c9SPeter Brune out from the matrix size which one it is. 825ab8d36c9SPeter Brune 826ab8d36c9SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 827ab8d36c9SPeter Brune is used. 828ab8d36c9SPeter Brune 829f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetInjection()` 830ab8d36c9SPeter Brune @*/ 831d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) 832d71ae5a4SJacob Faibussowitsch { 83322d28d08SBarry Smith SNES_FAS *fas; 834ab8d36c9SPeter Brune SNES levelsnes; 83522d28d08SBarry Smith 836ab8d36c9SPeter Brune PetscFunctionBegin; 837f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 838f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 8399566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 840ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 8419566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 8429566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->restrct)); 843ab8d36c9SPeter Brune fas->restrct = mat; 8443ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 845ab8d36c9SPeter Brune } 846ab8d36c9SPeter Brune 847ab8d36c9SPeter Brune /*@ 848ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the 849ab8d36c9SPeter Brune restriction from l to the l-1th level 850ab8d36c9SPeter Brune 851ab8d36c9SPeter Brune Input Parameters: 852f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 853ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 854ab8d36c9SPeter Brune 855f6dfbefdSBarry Smith Output Parameter: 856ab8d36c9SPeter Brune . mat - the interpolation operator 857ab8d36c9SPeter Brune 858ab8d36c9SPeter Brune Level: advanced 859ab8d36c9SPeter Brune 860f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetRestriction()`, `SNESFASGetInjection()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()` 861ab8d36c9SPeter Brune @*/ 862d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) 863d71ae5a4SJacob Faibussowitsch { 86422d28d08SBarry Smith SNES_FAS *fas; 865ab8d36c9SPeter Brune SNES levelsnes; 86622d28d08SBarry Smith 867ab8d36c9SPeter Brune PetscFunctionBegin; 868f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 8694f572ea9SToby Isaac PetscAssertPointer(mat, 3); 8709566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 871ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 872ab8d36c9SPeter Brune *mat = fas->restrct; 8733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 874ab8d36c9SPeter Brune } 875ab8d36c9SPeter Brune 876ab8d36c9SPeter Brune /*@ 877ab8d36c9SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 878ab8d36c9SPeter Brune from level l to l-1. 879ab8d36c9SPeter Brune 880ab8d36c9SPeter Brune Input Parameters: 881f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 882ab8d36c9SPeter Brune . mat - the restriction matrix 883ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 884ab8d36c9SPeter Brune 885ab8d36c9SPeter Brune Level: advanced 886ab8d36c9SPeter Brune 887f6dfbefdSBarry Smith Note: 888ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to 889ab8d36c9SPeter Brune project the solution instead. 890ab8d36c9SPeter Brune 891f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetRestriction()` 892ab8d36c9SPeter Brune @*/ 893d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) 894d71ae5a4SJacob Faibussowitsch { 89522d28d08SBarry Smith SNES_FAS *fas; 896ab8d36c9SPeter Brune SNES levelsnes; 89722d28d08SBarry Smith 898ab8d36c9SPeter Brune PetscFunctionBegin; 899f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 900f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 9019566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 902ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 9039566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 9049566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->inject)); 9051aa26658SKarl Rupp 906ab8d36c9SPeter Brune fas->inject = mat; 9073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 908ab8d36c9SPeter Brune } 909ab8d36c9SPeter Brune 910ab8d36c9SPeter Brune /*@ 911ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the 912ab8d36c9SPeter Brune injection from l-1 to the lth level 913ab8d36c9SPeter Brune 914ab8d36c9SPeter Brune Input Parameters: 915f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 916ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 917ab8d36c9SPeter Brune 918f6dfbefdSBarry Smith Output Parameter: 919ab8d36c9SPeter Brune . mat - the injection operator 920ab8d36c9SPeter Brune 921ab8d36c9SPeter Brune Level: advanced 922ab8d36c9SPeter Brune 923f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()` 924ab8d36c9SPeter Brune @*/ 925d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) 926d71ae5a4SJacob Faibussowitsch { 92722d28d08SBarry Smith SNES_FAS *fas; 928ab8d36c9SPeter Brune SNES levelsnes; 92922d28d08SBarry Smith 930ab8d36c9SPeter Brune PetscFunctionBegin; 931f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 9324f572ea9SToby Isaac PetscAssertPointer(mat, 3); 9339566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 934ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 935ab8d36c9SPeter Brune *mat = fas->inject; 9363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 937ab8d36c9SPeter Brune } 938ab8d36c9SPeter Brune 939ab8d36c9SPeter Brune /*@ 940ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 941ab8d36c9SPeter Brune operator from level l to l-1. 942ab8d36c9SPeter Brune 943ab8d36c9SPeter Brune Input Parameters: 944f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 945ab8d36c9SPeter Brune . rscale - the restriction scaling 946ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 947ab8d36c9SPeter Brune 948ab8d36c9SPeter Brune Level: advanced 949ab8d36c9SPeter Brune 950f6dfbefdSBarry Smith Note: 951ab8d36c9SPeter Brune This is only used in the case that the injection is not set. 952ab8d36c9SPeter Brune 953f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 954ab8d36c9SPeter Brune @*/ 955d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) 956d71ae5a4SJacob Faibussowitsch { 957ab8d36c9SPeter Brune SNES_FAS *fas; 958ab8d36c9SPeter Brune SNES levelsnes; 95922d28d08SBarry Smith 960ab8d36c9SPeter Brune PetscFunctionBegin; 961f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 962f833ba53SLisandro Dalcin if (rscale) PetscValidHeaderSpecific(rscale, VEC_CLASSID, 3); 9639566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 964ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 9659566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rscale)); 9669566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->rscale)); 967ab8d36c9SPeter Brune fas->rscale = rscale; 9683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 969ab8d36c9SPeter Brune } 970ab8d36c9SPeter Brune 971ab8d36c9SPeter Brune /*@ 972ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level. 973ab8d36c9SPeter Brune 974ab8d36c9SPeter Brune Input Parameters: 975f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 976ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 977ab8d36c9SPeter Brune 978f6dfbefdSBarry Smith Output Parameter: 979*ceaaa498SBarry Smith . smooth - the smoother 980ab8d36c9SPeter Brune 981ab8d36c9SPeter Brune Level: advanced 982ab8d36c9SPeter Brune 983f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 984ab8d36c9SPeter Brune @*/ 985d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) 986d71ae5a4SJacob Faibussowitsch { 987ab8d36c9SPeter Brune SNES_FAS *fas; 988ab8d36c9SPeter Brune SNES levelsnes; 98922d28d08SBarry Smith 990ab8d36c9SPeter Brune PetscFunctionBegin; 991f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 9924f572ea9SToby Isaac PetscAssertPointer(smooth, 3); 9939566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 994ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 99548a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 996ab8d36c9SPeter Brune *smooth = fas->smoothd; 9973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 998ab8d36c9SPeter Brune } 999ab8d36c9SPeter Brune 1000ab8d36c9SPeter Brune /*@ 1001ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level. 1002ab8d36c9SPeter Brune 1003ab8d36c9SPeter Brune Input Parameters: 1004f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1005ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1006ab8d36c9SPeter Brune 1007f6dfbefdSBarry Smith Output Parameter: 1008*ceaaa498SBarry Smith . smooth - the smoother 1009ab8d36c9SPeter Brune 1010ab8d36c9SPeter Brune Level: advanced 1011ab8d36c9SPeter Brune 1012f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1013ab8d36c9SPeter Brune @*/ 1014d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) 1015d71ae5a4SJacob Faibussowitsch { 1016ab8d36c9SPeter Brune SNES_FAS *fas; 1017ab8d36c9SPeter Brune SNES levelsnes; 101822d28d08SBarry Smith 1019ab8d36c9SPeter Brune PetscFunctionBegin; 1020f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 10214f572ea9SToby Isaac PetscAssertPointer(smooth, 3); 10229566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 1023ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1024ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 102548a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 102648a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu)); 1027ab8d36c9SPeter Brune *smooth = fas->smoothd; 10283ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1029ab8d36c9SPeter Brune } 1030ab8d36c9SPeter Brune 1031ab8d36c9SPeter Brune /*@ 1032ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level. 1033ab8d36c9SPeter Brune 1034ab8d36c9SPeter Brune Input Parameters: 1035f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1036ab8d36c9SPeter Brune - level - the level (0 is coarsest) 1037ab8d36c9SPeter Brune 1038f6dfbefdSBarry Smith Output Parameter: 1039*ceaaa498SBarry Smith . smooth - the smoother 1040ab8d36c9SPeter Brune 1041ab8d36c9SPeter Brune Level: advanced 1042ab8d36c9SPeter Brune 1043f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1044ab8d36c9SPeter Brune @*/ 1045d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) 1046d71ae5a4SJacob Faibussowitsch { 1047ab8d36c9SPeter Brune SNES_FAS *fas; 1048ab8d36c9SPeter Brune SNES levelsnes; 104922d28d08SBarry Smith 1050ab8d36c9SPeter Brune PetscFunctionBegin; 1051f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 10524f572ea9SToby Isaac PetscAssertPointer(smooth, 3); 10539566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 1054ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1055ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 105648a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 105748a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu)); 1058ab8d36c9SPeter Brune *smooth = fas->smoothu; 10593ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1060ab8d36c9SPeter Brune } 1061d6ad1212SPeter Brune 1062d6ad1212SPeter Brune /*@ 1063d6ad1212SPeter Brune SNESFASGetCoarseSolve - Gets the coarsest solver. 1064d6ad1212SPeter Brune 1065f6dfbefdSBarry Smith Input Parameter: 1066f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 1067d6ad1212SPeter Brune 1068f6dfbefdSBarry Smith Output Parameter: 1069a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver 1070d6ad1212SPeter Brune 1071d6ad1212SPeter Brune Level: advanced 1072d6ad1212SPeter Brune 1073f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1074d6ad1212SPeter Brune @*/ 1075d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse) 1076d71ae5a4SJacob Faibussowitsch { 1077d6ad1212SPeter Brune SNES_FAS *fas; 1078d6ad1212SPeter Brune SNES levelsnes; 107922d28d08SBarry Smith 1080d6ad1212SPeter Brune PetscFunctionBegin; 1081f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 10824f572ea9SToby Isaac PetscAssertPointer(coarse, 2); 10839566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, 0, &levelsnes)); 1084d6ad1212SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1085d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 108648a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 1087a3a80b83SMatthew G. Knepley *coarse = fas->smoothd; 10883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1089d6ad1212SPeter Brune } 1090928e959bSPeter Brune 1091928e959bSPeter Brune /*@ 1092f6dfbefdSBarry Smith SNESFASFullSetDownSweep - Smooth during the initial downsweep for `SNESFAS` 1093928e959bSPeter Brune 1094c3339decSBarry Smith Logically Collective 1095928e959bSPeter Brune 1096928e959bSPeter Brune Input Parameters: 1097f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1098928e959bSPeter Brune - swp - whether to downsweep or not 1099928e959bSPeter Brune 1100928e959bSPeter Brune Options Database Key: 1101928e959bSPeter Brune . -snes_fas_full_downsweep - Sets number of pre-smoothing steps 1102928e959bSPeter Brune 1103928e959bSPeter Brune Level: advanced 1104928e959bSPeter Brune 1105f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 1106928e959bSPeter Brune @*/ 1107d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullSetDownSweep(SNES snes, PetscBool swp) 1108d71ae5a4SJacob Faibussowitsch { 1109f833ba53SLisandro Dalcin SNES_FAS *fas; 1110928e959bSPeter Brune 1111928e959bSPeter Brune PetscFunctionBegin; 1112f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1113f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 1114928e959bSPeter Brune fas->full_downsweep = swp; 11151baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetDownSweep(fas->next, swp)); 11163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1117928e959bSPeter Brune } 11186dfbafefSToby Isaac 11196dfbafefSToby Isaac /*@ 11206dfbafefSToby Isaac SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 11216dfbafefSToby Isaac 1122c3339decSBarry Smith Logically Collective 11236dfbafefSToby Isaac 11246dfbafefSToby Isaac Input Parameters: 1125f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 11266dfbafefSToby Isaac - total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 11276dfbafefSToby Isaac 11286dfbafefSToby Isaac Options Database Key: 112967b8a455SSatish Balay . -snes_fas_full_total - Use total restriction and interpolation on the initial down and up sweeps for the full FAS cycle 11306dfbafefSToby Isaac 11316dfbafefSToby Isaac Level: advanced 11326dfbafefSToby Isaac 1133f6dfbefdSBarry Smith Note: 1134f6dfbefdSBarry Smith This option is only significant if the interpolation of a coarse correction (`MatInterpolate()`) is significantly different from total 1135f6dfbefdSBarry Smith solution interpolation (`DMInterpolateSolution()`). 11366dfbafefSToby Isaac 1137f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()` 11386dfbafefSToby Isaac @*/ 1139d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullSetTotal(SNES snes, PetscBool total) 1140d71ae5a4SJacob Faibussowitsch { 11416dfbafefSToby Isaac SNES_FAS *fas; 11426dfbafefSToby Isaac 11436dfbafefSToby Isaac PetscFunctionBegin; 11446dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 11456dfbafefSToby Isaac fas = (SNES_FAS *)snes->data; 11466dfbafefSToby Isaac fas->full_total = total; 11471baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetTotal(fas->next, total)); 11483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11496dfbafefSToby Isaac } 11506dfbafefSToby Isaac 11516dfbafefSToby Isaac /*@ 11526dfbafefSToby Isaac SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 11536dfbafefSToby Isaac 1154c3339decSBarry Smith Logically Collective 11556dfbafefSToby Isaac 1156f6dfbefdSBarry Smith Input Parameter: 1157f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 11586dfbafefSToby Isaac 1159*ceaaa498SBarry Smith Output Parameter: 11606dfbafefSToby Isaac . total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 11616dfbafefSToby Isaac 11626dfbafefSToby Isaac Level: advanced 11636dfbafefSToby Isaac 1164f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`, `SNESFullSetTotal()` 11656dfbafefSToby Isaac @*/ 1166d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullGetTotal(SNES snes, PetscBool *total) 1167d71ae5a4SJacob Faibussowitsch { 11686dfbafefSToby Isaac SNES_FAS *fas; 11696dfbafefSToby Isaac 11706dfbafefSToby Isaac PetscFunctionBegin; 11716dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 11726dfbafefSToby Isaac fas = (SNES_FAS *)snes->data; 11736dfbafefSToby Isaac *total = fas->full_total; 11743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11756dfbafefSToby Isaac } 1176