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); 50ab8d36c9SPeter Brune PetscValidPointer(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 132*e4094ef1SJacob 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); 145064a246eSJacob Faibussowitsch PetscValidIntPointer(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 /*@ 152f6dfbefdSBarry Smith SNESFASGetCycleSNES - Gets the `SNES` corresponding to a particular 153*e4094ef1SJacob Faibussowitsch 154*e4094ef1SJacob Faibussowitsch Level Of The `Snesfas` Hierarchy. 155ab8d36c9SPeter Brune 156ab8d36c9SPeter Brune Input Parameters: 157f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context 158f6dfbefdSBarry Smith - level - the level to get 159f6dfbefdSBarry Smith 160f6dfbefdSBarry Smith Output Parameter: 161f6dfbefdSBarry Smith . lsnes - the `SNES` for the requested level 162ab8d36c9SPeter Brune 163ab8d36c9SPeter Brune Level: advanced 164ab8d36c9SPeter Brune 165f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()`, `SNESFASGetLevels()` 166ab8d36c9SPeter Brune @*/ 167d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetCycleSNES(SNES snes, PetscInt level, SNES *lsnes) 168d71ae5a4SJacob Faibussowitsch { 169f833ba53SLisandro Dalcin SNES_FAS *fas; 170ab8d36c9SPeter Brune PetscInt i; 171ab8d36c9SPeter Brune 172ab8d36c9SPeter Brune PetscFunctionBegin; 173f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 174f833ba53SLisandro Dalcin PetscValidPointer(lsnes, 3); 175f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 17663a3b9bcSJacob 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); 1770b121fc5SBarry Smith PetscCheck(fas->level == fas->levels - 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetCycleSNES may only be called on the finest-level SNES"); 178ab8d36c9SPeter Brune 179ab8d36c9SPeter Brune *lsnes = snes; 180ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) { 181ab8d36c9SPeter Brune *lsnes = fas->next; 182ab8d36c9SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 183ab8d36c9SPeter Brune } 18408401ef6SPierre Jolivet PetscCheck(fas->level == level, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "SNESFAS level hierarchy corrupt"); 1853ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 186ab8d36c9SPeter Brune } 187ab8d36c9SPeter Brune 188ab8d36c9SPeter Brune /*@ 189ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 190ab8d36c9SPeter Brune use on all levels. 191ab8d36c9SPeter Brune 192c3339decSBarry Smith Logically Collective 193ab8d36c9SPeter Brune 194ab8d36c9SPeter Brune Input Parameters: 195f6dfbefdSBarry Smith + snes - the `SNES` nonlinear multigrid context 196ab8d36c9SPeter Brune - n - the number of smoothing steps 197ab8d36c9SPeter Brune 198ab8d36c9SPeter Brune Options Database Key: 199ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 200ab8d36c9SPeter Brune 201ab8d36c9SPeter Brune Level: advanced 202ab8d36c9SPeter Brune 203f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothDown()` 204ab8d36c9SPeter Brune @*/ 205d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) 206d71ae5a4SJacob Faibussowitsch { 207f833ba53SLisandro Dalcin SNES_FAS *fas; 20822d28d08SBarry Smith 209ab8d36c9SPeter Brune PetscFunctionBegin; 210f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 211f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 212ab8d36c9SPeter Brune fas->max_up_it = n; 21348a46eb9SPierre Jolivet if (!fas->smoothu && fas->level != 0) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu)); 2141baa6e33SBarry Smith if (fas->smoothu) PetscCall(SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs)); 2151baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothUp(fas->next, n)); 2163ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 217ab8d36c9SPeter Brune } 218ab8d36c9SPeter Brune 219ab8d36c9SPeter Brune /*@ 220ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 221ab8d36c9SPeter Brune use on all levels. 222ab8d36c9SPeter Brune 223c3339decSBarry Smith Logically Collective 224ab8d36c9SPeter Brune 225ab8d36c9SPeter Brune Input Parameters: 226f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 227ab8d36c9SPeter Brune - n - the number of smoothing steps 228ab8d36c9SPeter Brune 229ab8d36c9SPeter Brune Options Database Key: 230ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 231ab8d36c9SPeter Brune 232ab8d36c9SPeter Brune Level: advanced 233ab8d36c9SPeter Brune 234f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 235ab8d36c9SPeter Brune @*/ 236d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) 237d71ae5a4SJacob Faibussowitsch { 238f833ba53SLisandro Dalcin SNES_FAS *fas; 23922d28d08SBarry Smith 240ab8d36c9SPeter Brune PetscFunctionBegin; 241f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 242f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 24348a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd)); 2449566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs)); 2451aa26658SKarl Rupp 246ab8d36c9SPeter Brune fas->max_down_it = n; 2471baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetNumberSmoothDown(fas->next, n)); 2483ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 249ab8d36c9SPeter Brune } 250ab8d36c9SPeter Brune 25187f44e3fSPeter Brune /*@ 252f6dfbefdSBarry Smith SNESFASSetContinuation - Sets the `SNESFAS` cycle to default to exact Newton solves on the upsweep 25387f44e3fSPeter Brune 254c3339decSBarry Smith Logically Collective 25587f44e3fSPeter Brune 25687f44e3fSPeter Brune Input Parameters: 257f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 258*e4094ef1SJacob Faibussowitsch - continuation - the number of smoothing steps 25987f44e3fSPeter Brune 26087f44e3fSPeter Brune Options Database Key: 26187f44e3fSPeter Brune . -snes_fas_continuation - sets continuation to true 26287f44e3fSPeter Brune 26387f44e3fSPeter Brune Level: advanced 26487f44e3fSPeter Brune 265f6dfbefdSBarry Smith Note: 26695452b02SPatrick Sanan This sets the prefix on the upsweep smoothers to -fas_continuation 26787f44e3fSPeter Brune 268f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 26987f44e3fSPeter Brune @*/ 270d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetContinuation(SNES snes, PetscBool continuation) 271d71ae5a4SJacob Faibussowitsch { 27287f44e3fSPeter Brune const char *optionsprefix; 27387f44e3fSPeter Brune char tprefix[128]; 274f833ba53SLisandro Dalcin SNES_FAS *fas; 27587f44e3fSPeter Brune 27687f44e3fSPeter Brune PetscFunctionBegin; 277f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 278f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 2799566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 28048a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu)); 2819566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_levels_continuation_", sizeof(tprefix))); 2829566063dSJacob Faibussowitsch PetscCall(SNESSetOptionsPrefix(fas->smoothu, optionsprefix)); 2839566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(fas->smoothu, tprefix)); 2849566063dSJacob Faibussowitsch PetscCall(SNESSetType(fas->smoothu, SNESNEWTONLS)); 2859566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(fas->smoothu, fas->fine->abstol, fas->fine->rtol, fas->fine->stol, 50, 100)); 28687f44e3fSPeter Brune fas->continuation = continuation; 2871baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetContinuation(fas->next, continuation)); 2883ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 28987f44e3fSPeter Brune } 29087f44e3fSPeter Brune 291ab8d36c9SPeter Brune /*@ 292f6dfbefdSBarry Smith SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use `SNESFASSetCyclesOnLevel()` for more 293ab8d36c9SPeter Brune complicated cycling. 294ab8d36c9SPeter Brune 295c3339decSBarry Smith Logically Collective 296ab8d36c9SPeter Brune 297ab8d36c9SPeter Brune Input Parameters: 298f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 299ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 300ab8d36c9SPeter Brune 301ab8d36c9SPeter Brune Options Database Key: 30267b8a455SSatish Balay . -snes_fas_cycles <1,2> - 1 for V-cycle, 2 for W-cycle 303ab8d36c9SPeter Brune 304ab8d36c9SPeter Brune Level: advanced 305ab8d36c9SPeter Brune 306f6dfbefdSBarry Smith .seealso: `SNES`, `SNESFAS`, `SNESFASSetCyclesOnLevel()` 307ab8d36c9SPeter Brune @*/ 308d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) 309d71ae5a4SJacob Faibussowitsch { 310f833ba53SLisandro Dalcin SNES_FAS *fas; 311ab8d36c9SPeter Brune PetscBool isFine; 31222d28d08SBarry Smith 313ab8d36c9SPeter Brune PetscFunctionBegin; 314f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3159566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 316f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 317ab8d36c9SPeter Brune fas->n_cycles = cycles; 31848a46eb9SPierre Jolivet if (!isFine) PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs)); 3191baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASSetCycles(fas->next, cycles)); 3203ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 321ab8d36c9SPeter Brune } 322ab8d36c9SPeter Brune 323c8c899caSPeter Brune /*@ 324c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring 325c8c899caSPeter Brune 326c3339decSBarry Smith Logically Collective 327c8c899caSPeter Brune 328c8c899caSPeter Brune Input Parameters: 329f6dfbefdSBarry Smith + snes - the `SNESFAS` context 3302fe279fdSBarry Smith . vf - viewer and format structure (may be `NULL` if flg is `PETSC_FALSE`) 331c8c899caSPeter Brune - flg - monitor or not 332c8c899caSPeter Brune 333c8c899caSPeter Brune Level: advanced 334c8c899caSPeter Brune 335f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESSetMonitor()`, `SNESFASSetCyclesOnLevel()` 336c8c899caSPeter Brune @*/ 337d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg) 338d71ae5a4SJacob Faibussowitsch { 339f833ba53SLisandro Dalcin SNES_FAS *fas; 340c8c899caSPeter Brune PetscBool isFine; 341f833ba53SLisandro Dalcin PetscInt i, levels; 342c8c899caSPeter Brune SNES levelsnes; 34322d28d08SBarry Smith 344c8c899caSPeter Brune PetscFunctionBegin; 345f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3469566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 347f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 348f833ba53SLisandro Dalcin levels = fas->levels; 349c8c899caSPeter Brune if (isFine) { 350c8c899caSPeter Brune for (i = 0; i < levels; i++) { 3519566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 352c8c899caSPeter Brune fas = (SNES_FAS *)levelsnes->data; 353c8c899caSPeter Brune if (flg) { 354c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */ 3559566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes)); 356d142ab34SLawrence Mitchell /* Only register destroy on finest level */ 3579566063dSJacob Faibussowitsch PetscCall(SNESMonitorSet(levelsnes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))SNESMonitorDefault, vf, (!i ? (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy : NULL))); 3581aa26658SKarl Rupp } else if (i != fas->levels - 1) { 359c8c899caSPeter Brune /* unset the monitors on the coarse levels */ 3609566063dSJacob Faibussowitsch PetscCall(SNESMonitorCancel(levelsnes)); 361c8c899caSPeter Brune } 362c8c899caSPeter Brune } 363c8c899caSPeter Brune } 3643ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 365c8c899caSPeter Brune } 366c8c899caSPeter Brune 3670dd27c6cSPeter Brune /*@ 368f6dfbefdSBarry Smith SNESFASSetLog - Sets or unsets time logging for various `SNESFAS` stages on all levels 3690dd27c6cSPeter Brune 370c3339decSBarry Smith Logically Collective 3710dd27c6cSPeter Brune 3720dd27c6cSPeter Brune Input Parameters: 373f6dfbefdSBarry Smith + snes - the `SNESFAS` context 3740dd27c6cSPeter Brune - flg - monitor or not 3750dd27c6cSPeter Brune 3760dd27c6cSPeter Brune Level: advanced 3770dd27c6cSPeter Brune 378f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetMonitor()` 3790dd27c6cSPeter Brune @*/ 380d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) 381d71ae5a4SJacob Faibussowitsch { 382f833ba53SLisandro Dalcin SNES_FAS *fas; 3830dd27c6cSPeter Brune PetscBool isFine; 384f833ba53SLisandro Dalcin PetscInt i, levels; 3850dd27c6cSPeter Brune SNES levelsnes; 3860dd27c6cSPeter Brune char eventname[128]; 3870dd27c6cSPeter Brune 3880dd27c6cSPeter Brune PetscFunctionBegin; 389f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 3909566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 391f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 392f833ba53SLisandro Dalcin levels = fas->levels; 3930dd27c6cSPeter Brune if (isFine) { 3940dd27c6cSPeter Brune for (i = 0; i < levels; i++) { 3959566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 3960dd27c6cSPeter Brune fas = (SNES_FAS *)levelsnes->data; 3970dd27c6cSPeter Brune if (flg) { 3989566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSetup %d", (int)i)); 3999566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsetup)); 4009566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSmooth %d", (int)i)); 4019566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsolve)); 4029566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASResid %d", (int)i)); 4039566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventresidual)); 4049566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASInterp %d", (int)i)); 4059566063dSJacob Faibussowitsch PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventinterprestrict)); 4060dd27c6cSPeter Brune } else { 4070298fd71SBarry Smith fas->eventsmoothsetup = 0; 4080298fd71SBarry Smith fas->eventsmoothsolve = 0; 4090298fd71SBarry Smith fas->eventresidual = 0; 4100298fd71SBarry Smith fas->eventinterprestrict = 0; 4110dd27c6cSPeter Brune } 4120dd27c6cSPeter Brune } 4130dd27c6cSPeter Brune } 4143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 4150dd27c6cSPeter Brune } 4160dd27c6cSPeter Brune 417ab8d36c9SPeter Brune /* 418ab8d36c9SPeter Brune Creates the default smoother type. 419ab8d36c9SPeter Brune 42004d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level. 421ab8d36c9SPeter Brune 422ab8d36c9SPeter Brune */ 423d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) 424d71ae5a4SJacob Faibussowitsch { 425ab8d36c9SPeter Brune SNES_FAS *fas; 426ab8d36c9SPeter Brune const char *optionsprefix; 427ab8d36c9SPeter Brune char tprefix[128]; 428ab8d36c9SPeter Brune SNES nsmooth; 42922d28d08SBarry Smith 430ab8d36c9SPeter Brune PetscFunctionBegin; 431f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 432064a246eSJacob Faibussowitsch PetscValidPointer(smooth, 2); 433ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 4349566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix)); 435ab8d36c9SPeter Brune /* create the default smoother */ 4369566063dSJacob Faibussowitsch PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth)); 437ab8d36c9SPeter Brune if (fas->level == 0) { 4389566063dSJacob Faibussowitsch PetscCall(PetscStrncpy(tprefix, "fas_coarse_", sizeof(tprefix))); 4399566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix)); 4409566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix)); 4419566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNEWTONLS)); 4429566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs)); 443ab8d36c9SPeter Brune } else { 4449566063dSJacob Faibussowitsch PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%d_", (int)fas->level)); 4459566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix)); 4469566063dSJacob Faibussowitsch PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix)); 4479566063dSJacob Faibussowitsch PetscCall(SNESSetType(nsmooth, SNESNRICHARDSON)); 4489566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs)); 449ab8d36c9SPeter Brune } 4509566063dSJacob Faibussowitsch PetscCall(PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1)); 4519566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth)); 4529566063dSJacob Faibussowitsch PetscCall(PetscObjectComposedDataSetInt((PetscObject)nsmooth, PetscMGLevelId, fas->level)); 453ab8d36c9SPeter Brune *smooth = nsmooth; 4543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 455ab8d36c9SPeter Brune } 456ab8d36c9SPeter Brune 457ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */ 458ab8d36c9SPeter Brune 459ab8d36c9SPeter Brune /*@ 460ab8d36c9SPeter Brune SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 461ab8d36c9SPeter Brune 462c3339decSBarry Smith Logically Collective 463ab8d36c9SPeter Brune 464ab8d36c9SPeter Brune Input Parameters: 465f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 466ab8d36c9SPeter Brune . level - the level to set the number of cycles on 467ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 468ab8d36c9SPeter Brune 469ab8d36c9SPeter Brune Level: advanced 470ab8d36c9SPeter Brune 471f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetCycles()` 472ab8d36c9SPeter Brune @*/ 473d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) 474d71ae5a4SJacob Faibussowitsch { 475f833ba53SLisandro Dalcin SNES_FAS *fas; 47622d28d08SBarry Smith 477ab8d36c9SPeter Brune PetscFunctionBegin; 478f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 479f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 480ab8d36c9SPeter Brune fas->n_cycles = cycles; 4819566063dSJacob Faibussowitsch PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs)); 4823ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 483ab8d36c9SPeter Brune } 484ab8d36c9SPeter Brune 485ab8d36c9SPeter Brune /*@ 486ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 487ab8d36c9SPeter Brune 488c3339decSBarry Smith Logically Collective 489ab8d36c9SPeter Brune 490f6dfbefdSBarry Smith Input Parameter: 491f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 492ab8d36c9SPeter Brune 493f6dfbefdSBarry Smith Output Parameter: 494ab8d36c9SPeter Brune . smooth - the smoother 495ab8d36c9SPeter Brune 496ab8d36c9SPeter Brune Level: advanced 497ab8d36c9SPeter Brune 498f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()` 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); 506f833ba53SLisandro Dalcin PetscValidPointer(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 522f6dfbefdSBarry Smith Note: 523ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent 524ab8d36c9SPeter Brune default behavior in the process of the solve. 525ab8d36c9SPeter Brune 526ab8d36c9SPeter Brune Level: advanced 527ab8d36c9SPeter Brune 528f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherDown()` 529ab8d36c9SPeter Brune @*/ 530d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 531d71ae5a4SJacob Faibussowitsch { 532ab8d36c9SPeter Brune SNES_FAS *fas; 5335fd66863SKarl Rupp 534ab8d36c9SPeter Brune PetscFunctionBegin; 535f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 536f833ba53SLisandro Dalcin PetscValidPointer(smoothu, 2); 537ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 5381aa26658SKarl Rupp if (!fas->smoothu) *smoothu = fas->smoothd; 5391aa26658SKarl Rupp else *smoothu = fas->smoothu; 5403ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 541ab8d36c9SPeter Brune } 542ab8d36c9SPeter Brune 543ab8d36c9SPeter Brune /*@ 544ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 545ab8d36c9SPeter Brune 546c3339decSBarry Smith Logically Collective 547ab8d36c9SPeter Brune 548f6dfbefdSBarry Smith Input Parameter: 549f6dfbefdSBarry Smith . snes - `SNESFAS`, the nonlinear multigrid context 550ab8d36c9SPeter Brune 551f6dfbefdSBarry Smith Output Parameter: 552ab8d36c9SPeter Brune . smoothd - the smoother 553ab8d36c9SPeter Brune 554ab8d36c9SPeter Brune Level: advanced 555ab8d36c9SPeter Brune 556f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 557ab8d36c9SPeter Brune @*/ 558d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 559d71ae5a4SJacob Faibussowitsch { 560ab8d36c9SPeter Brune SNES_FAS *fas; 5615fd66863SKarl Rupp 562ab8d36c9SPeter Brune PetscFunctionBegin; 563f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 564f833ba53SLisandro Dalcin PetscValidPointer(smoothd, 2); 565ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 566ab8d36c9SPeter Brune *smoothd = fas->smoothd; 5673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 568ab8d36c9SPeter Brune } 569ab8d36c9SPeter Brune 570ab8d36c9SPeter Brune /*@ 571ab8d36c9SPeter Brune SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 572ab8d36c9SPeter Brune 573c3339decSBarry Smith Logically Collective 574ab8d36c9SPeter Brune 575f6dfbefdSBarry Smith Input Parameter: 576f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 577ab8d36c9SPeter Brune 578f6dfbefdSBarry Smith Output Parameter: 579f6dfbefdSBarry Smith . correction - the coarse correction solve on this level 580ab8d36c9SPeter Brune 581f6dfbefdSBarry Smith Note: 5820298fd71SBarry Smith Returns NULL on the coarsest level. 583ab8d36c9SPeter Brune 584ab8d36c9SPeter Brune Level: advanced 585ab8d36c9SPeter Brune 586f6dfbefdSBarry Smith .seealso: `SNESFAS` `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 587ab8d36c9SPeter Brune @*/ 588d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 589d71ae5a4SJacob Faibussowitsch { 590ab8d36c9SPeter Brune SNES_FAS *fas; 5915fd66863SKarl Rupp 592ab8d36c9SPeter Brune PetscFunctionBegin; 593f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 594f833ba53SLisandro Dalcin PetscValidPointer(correction, 2); 595ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 596ab8d36c9SPeter Brune *correction = fas->next; 5973ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 598ab8d36c9SPeter Brune } 599ab8d36c9SPeter Brune 600ab8d36c9SPeter Brune /*@ 601ab8d36c9SPeter Brune SNESFASCycleGetInterpolation - Gets the interpolation on this level 602ab8d36c9SPeter Brune 603c3339decSBarry Smith Logically Collective 604ab8d36c9SPeter Brune 605f6dfbefdSBarry Smith Input Parameter: 606f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 607ab8d36c9SPeter Brune 608f6dfbefdSBarry Smith Output Parameter: 609ab8d36c9SPeter Brune . mat - the interpolation operator on this level 610ab8d36c9SPeter Brune 611f6dfbefdSBarry Smith Level: advanced 612ab8d36c9SPeter Brune 613f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()` 614ab8d36c9SPeter Brune @*/ 615d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 616d71ae5a4SJacob Faibussowitsch { 617ab8d36c9SPeter Brune SNES_FAS *fas; 6185fd66863SKarl Rupp 619ab8d36c9SPeter Brune PetscFunctionBegin; 620f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 621f833ba53SLisandro Dalcin PetscValidPointer(mat, 2); 622ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 623ab8d36c9SPeter Brune *mat = fas->interpolate; 6243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 625ab8d36c9SPeter Brune } 626ab8d36c9SPeter Brune 627ab8d36c9SPeter Brune /*@ 628ab8d36c9SPeter Brune SNESFASCycleGetRestriction - Gets the restriction on this level 629ab8d36c9SPeter Brune 630c3339decSBarry Smith Logically Collective 631ab8d36c9SPeter Brune 632f6dfbefdSBarry Smith Input Parameter: 633f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 634ab8d36c9SPeter Brune 635f6dfbefdSBarry Smith Output Parameter: 636ab8d36c9SPeter Brune . mat - the restriction operator on this level 637ab8d36c9SPeter Brune 638f6dfbefdSBarry Smith Level: advanced 639ab8d36c9SPeter Brune 640f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetRestriction()`, `SNESFASCycleGetInterpolation()` 641ab8d36c9SPeter Brune @*/ 642d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 643d71ae5a4SJacob Faibussowitsch { 644ab8d36c9SPeter Brune SNES_FAS *fas; 6455fd66863SKarl Rupp 646ab8d36c9SPeter Brune PetscFunctionBegin; 647f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 648f833ba53SLisandro Dalcin PetscValidPointer(mat, 2); 649ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 650ab8d36c9SPeter Brune *mat = fas->restrct; 6513ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 652ab8d36c9SPeter Brune } 653ab8d36c9SPeter Brune 654ab8d36c9SPeter Brune /*@ 655ab8d36c9SPeter Brune SNESFASCycleGetInjection - Gets the injection on this level 656ab8d36c9SPeter Brune 657c3339decSBarry Smith Logically Collective 658ab8d36c9SPeter Brune 659f6dfbefdSBarry Smith Input Parameter: 660f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 661ab8d36c9SPeter Brune 662f6dfbefdSBarry Smith Output Parameter: 663ab8d36c9SPeter Brune . mat - the restriction operator on this level 664ab8d36c9SPeter Brune 665f6dfbefdSBarry Smith Level: advanced 666ab8d36c9SPeter Brune 667f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASGetInjection()`, `SNESFASCycleGetRestriction()` 668ab8d36c9SPeter Brune @*/ 669d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 670d71ae5a4SJacob Faibussowitsch { 671ab8d36c9SPeter Brune SNES_FAS *fas; 6725fd66863SKarl Rupp 673ab8d36c9SPeter Brune PetscFunctionBegin; 674f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 675f833ba53SLisandro Dalcin PetscValidPointer(mat, 2); 676ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 677ab8d36c9SPeter Brune *mat = fas->inject; 6783ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 679ab8d36c9SPeter Brune } 680ab8d36c9SPeter Brune 681ab8d36c9SPeter Brune /*@ 682ab8d36c9SPeter Brune SNESFASCycleGetRScale - Gets the injection on this level 683ab8d36c9SPeter Brune 684c3339decSBarry Smith Logically Collective 685ab8d36c9SPeter Brune 686f6dfbefdSBarry Smith Input Parameter: 687f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 688ab8d36c9SPeter Brune 689f6dfbefdSBarry Smith Output Parameter: 690*e4094ef1SJacob Faibussowitsch . vec - the restriction operator on this level 691ab8d36c9SPeter Brune 692f6dfbefdSBarry Smith Level: advanced 693ab8d36c9SPeter Brune 694f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASCycleGetRestriction()`, `SNESFASGetRScale()` 695ab8d36c9SPeter Brune @*/ 696d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 697d71ae5a4SJacob Faibussowitsch { 698ab8d36c9SPeter Brune SNES_FAS *fas; 6995fd66863SKarl Rupp 700ab8d36c9SPeter Brune PetscFunctionBegin; 701f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 702f833ba53SLisandro Dalcin PetscValidPointer(vec, 2); 703ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 704ab8d36c9SPeter Brune *vec = fas->rscale; 7053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 706ab8d36c9SPeter Brune } 707ab8d36c9SPeter Brune 708ab8d36c9SPeter Brune /*@ 709ab8d36c9SPeter Brune SNESFASCycleIsFine - Determines if a given cycle is the fine level. 710ab8d36c9SPeter Brune 711c3339decSBarry Smith Logically Collective 712ab8d36c9SPeter Brune 713f6dfbefdSBarry Smith Input Parameter: 714f6dfbefdSBarry Smith . snes - the `SNESFAS` `SNES` context 715ab8d36c9SPeter Brune 716f6dfbefdSBarry Smith Output Parameter: 717ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not 718ab8d36c9SPeter Brune 719ab8d36c9SPeter Brune Level: advanced 720ab8d36c9SPeter Brune 721f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetLevels()` 722ab8d36c9SPeter Brune @*/ 723d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 724d71ae5a4SJacob Faibussowitsch { 725ab8d36c9SPeter Brune SNES_FAS *fas; 7265fd66863SKarl Rupp 727ab8d36c9SPeter Brune PetscFunctionBegin; 728f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 729534a8f05SLisandro Dalcin PetscValidBoolPointer(flg, 2); 730ab8d36c9SPeter Brune fas = (SNES_FAS *)snes->data; 7311aa26658SKarl Rupp if (fas->level == fas->levels - 1) *flg = PETSC_TRUE; 7321aa26658SKarl Rupp else *flg = PETSC_FALSE; 7333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 734ab8d36c9SPeter Brune } 735ab8d36c9SPeter Brune 736f6dfbefdSBarry Smith /* functions called on the finest level that return level-specific information */ 737ab8d36c9SPeter Brune 738ab8d36c9SPeter Brune /*@ 739f6dfbefdSBarry Smith SNESFASSetInterpolation - Sets the `Mat` to be used to apply the 740ab8d36c9SPeter Brune interpolation from l-1 to the lth level 741ab8d36c9SPeter Brune 742ab8d36c9SPeter Brune Input Parameters: 743f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 744ab8d36c9SPeter Brune . mat - the interpolation operator 745ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 746ab8d36c9SPeter Brune 747ab8d36c9SPeter Brune Level: advanced 748ab8d36c9SPeter Brune 749ab8d36c9SPeter Brune Notes: 750ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction 751ab8d36c9SPeter Brune for the same level. 752ab8d36c9SPeter Brune 753ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 754ab8d36c9SPeter Brune out from the matrix size which one it is. 755ab8d36c9SPeter Brune 756f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()`, `SNESFASSetRScale()` 757ab8d36c9SPeter Brune @*/ 758d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) 759d71ae5a4SJacob Faibussowitsch { 76022d28d08SBarry Smith SNES_FAS *fas; 761ab8d36c9SPeter Brune SNES levelsnes; 76222d28d08SBarry Smith 763ab8d36c9SPeter Brune PetscFunctionBegin; 764f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 765f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 7669566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 767ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 7689566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 7699566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->interpolate)); 770ab8d36c9SPeter Brune fas->interpolate = mat; 7713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 772ab8d36c9SPeter Brune } 773ab8d36c9SPeter Brune 774ab8d36c9SPeter Brune /*@ 775ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the 776ab8d36c9SPeter Brune interpolation from l-1 to the lth level 777ab8d36c9SPeter Brune 778ab8d36c9SPeter Brune Input Parameters: 779f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 780ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 781ab8d36c9SPeter Brune 782f6dfbefdSBarry Smith Output Parameter: 783ab8d36c9SPeter Brune . mat - the interpolation operator 784ab8d36c9SPeter Brune 785ab8d36c9SPeter Brune Level: advanced 786ab8d36c9SPeter Brune 787f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASGetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetRScale()` 788ab8d36c9SPeter Brune @*/ 789d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) 790d71ae5a4SJacob Faibussowitsch { 79122d28d08SBarry Smith SNES_FAS *fas; 792ab8d36c9SPeter Brune SNES levelsnes; 79322d28d08SBarry Smith 794ab8d36c9SPeter Brune PetscFunctionBegin; 795f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 796f833ba53SLisandro Dalcin PetscValidPointer(mat, 3); 7979566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 798ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 799ab8d36c9SPeter Brune *mat = fas->interpolate; 8003ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 801ab8d36c9SPeter Brune } 802ab8d36c9SPeter Brune 803ab8d36c9SPeter Brune /*@ 804f6dfbefdSBarry Smith SNESFASSetRestriction - Sets the matrix to be used to restrict the defect 805ab8d36c9SPeter Brune from level l to l-1. 806ab8d36c9SPeter Brune 807ab8d36c9SPeter Brune Input Parameters: 808f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 809ab8d36c9SPeter Brune . mat - the restriction matrix 810ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 811ab8d36c9SPeter Brune 812ab8d36c9SPeter Brune Level: advanced 813ab8d36c9SPeter Brune 814ab8d36c9SPeter Brune Notes: 815ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation 816ab8d36c9SPeter Brune for the same level. 817ab8d36c9SPeter Brune 818ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 819ab8d36c9SPeter Brune out from the matrix size which one it is. 820ab8d36c9SPeter Brune 821ab8d36c9SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 822ab8d36c9SPeter Brune is used. 823ab8d36c9SPeter Brune 824f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetInjection()` 825ab8d36c9SPeter Brune @*/ 826d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) 827d71ae5a4SJacob Faibussowitsch { 82822d28d08SBarry Smith SNES_FAS *fas; 829ab8d36c9SPeter Brune SNES levelsnes; 83022d28d08SBarry Smith 831ab8d36c9SPeter Brune PetscFunctionBegin; 832f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 833f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 8349566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 835ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 8369566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 8379566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->restrct)); 838ab8d36c9SPeter Brune fas->restrct = mat; 8393ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 840ab8d36c9SPeter Brune } 841ab8d36c9SPeter Brune 842ab8d36c9SPeter Brune /*@ 843ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the 844ab8d36c9SPeter Brune restriction from l to the l-1th level 845ab8d36c9SPeter Brune 846ab8d36c9SPeter Brune Input Parameters: 847f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 848ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 849ab8d36c9SPeter Brune 850f6dfbefdSBarry Smith Output Parameter: 851ab8d36c9SPeter Brune . mat - the interpolation operator 852ab8d36c9SPeter Brune 853ab8d36c9SPeter Brune Level: advanced 854ab8d36c9SPeter Brune 855f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetRestriction()`, `SNESFASGetInjection()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()` 856ab8d36c9SPeter Brune @*/ 857d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) 858d71ae5a4SJacob Faibussowitsch { 85922d28d08SBarry Smith SNES_FAS *fas; 860ab8d36c9SPeter Brune SNES levelsnes; 86122d28d08SBarry Smith 862ab8d36c9SPeter Brune PetscFunctionBegin; 863f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 864f833ba53SLisandro Dalcin PetscValidPointer(mat, 3); 8659566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 866ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 867ab8d36c9SPeter Brune *mat = fas->restrct; 8683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 869ab8d36c9SPeter Brune } 870ab8d36c9SPeter Brune 871ab8d36c9SPeter Brune /*@ 872ab8d36c9SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 873ab8d36c9SPeter Brune from level l to l-1. 874ab8d36c9SPeter Brune 875ab8d36c9SPeter Brune Input Parameters: 876f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 877ab8d36c9SPeter Brune . mat - the restriction matrix 878ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 879ab8d36c9SPeter Brune 880ab8d36c9SPeter Brune Level: advanced 881ab8d36c9SPeter Brune 882f6dfbefdSBarry Smith Note: 883ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to 884ab8d36c9SPeter Brune project the solution instead. 885ab8d36c9SPeter Brune 886f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInterpolation()`, `SNESFASSetRestriction()` 887ab8d36c9SPeter Brune @*/ 888d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) 889d71ae5a4SJacob Faibussowitsch { 89022d28d08SBarry Smith SNES_FAS *fas; 891ab8d36c9SPeter Brune SNES levelsnes; 89222d28d08SBarry Smith 893ab8d36c9SPeter Brune PetscFunctionBegin; 894f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 895f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3); 8969566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 897ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 8989566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)mat)); 8999566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->inject)); 9001aa26658SKarl Rupp 901ab8d36c9SPeter Brune fas->inject = mat; 9023ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 903ab8d36c9SPeter Brune } 904ab8d36c9SPeter Brune 905ab8d36c9SPeter Brune /*@ 906ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the 907ab8d36c9SPeter Brune injection from l-1 to the lth level 908ab8d36c9SPeter Brune 909ab8d36c9SPeter Brune Input Parameters: 910f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 911ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 912ab8d36c9SPeter Brune 913f6dfbefdSBarry Smith Output Parameter: 914ab8d36c9SPeter Brune . mat - the injection operator 915ab8d36c9SPeter Brune 916ab8d36c9SPeter Brune Level: advanced 917ab8d36c9SPeter Brune 918f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()` 919ab8d36c9SPeter Brune @*/ 920d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) 921d71ae5a4SJacob Faibussowitsch { 92222d28d08SBarry Smith SNES_FAS *fas; 923ab8d36c9SPeter Brune SNES levelsnes; 92422d28d08SBarry Smith 925ab8d36c9SPeter Brune PetscFunctionBegin; 926f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 927f833ba53SLisandro Dalcin PetscValidPointer(mat, 3); 9289566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 929ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 930ab8d36c9SPeter Brune *mat = fas->inject; 9313ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 932ab8d36c9SPeter Brune } 933ab8d36c9SPeter Brune 934ab8d36c9SPeter Brune /*@ 935ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 936ab8d36c9SPeter Brune operator from level l to l-1. 937ab8d36c9SPeter Brune 938ab8d36c9SPeter Brune Input Parameters: 939f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 940ab8d36c9SPeter Brune . rscale - the restriction scaling 941ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 942ab8d36c9SPeter Brune 943ab8d36c9SPeter Brune Level: advanced 944ab8d36c9SPeter Brune 945f6dfbefdSBarry Smith Note: 946ab8d36c9SPeter Brune This is only used in the case that the injection is not set. 947ab8d36c9SPeter Brune 948f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 949ab8d36c9SPeter Brune @*/ 950d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) 951d71ae5a4SJacob Faibussowitsch { 952ab8d36c9SPeter Brune SNES_FAS *fas; 953ab8d36c9SPeter Brune SNES levelsnes; 95422d28d08SBarry Smith 955ab8d36c9SPeter Brune PetscFunctionBegin; 956f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 957f833ba53SLisandro Dalcin if (rscale) PetscValidHeaderSpecific(rscale, VEC_CLASSID, 3); 9589566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 959ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 9609566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)rscale)); 9619566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->rscale)); 962ab8d36c9SPeter Brune fas->rscale = rscale; 9633ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 964ab8d36c9SPeter Brune } 965ab8d36c9SPeter Brune 966ab8d36c9SPeter Brune /*@ 967ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level. 968ab8d36c9SPeter Brune 969ab8d36c9SPeter Brune Input Parameters: 970f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 971ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 972ab8d36c9SPeter Brune 973f6dfbefdSBarry Smith Output Parameter: 974ab8d36c9SPeter Brune smooth - the smoother 975ab8d36c9SPeter Brune 976ab8d36c9SPeter Brune Level: advanced 977ab8d36c9SPeter Brune 978f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 979ab8d36c9SPeter Brune @*/ 980d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) 981d71ae5a4SJacob Faibussowitsch { 982ab8d36c9SPeter Brune SNES_FAS *fas; 983ab8d36c9SPeter Brune SNES levelsnes; 98422d28d08SBarry Smith 985ab8d36c9SPeter Brune PetscFunctionBegin; 986f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 987f833ba53SLisandro Dalcin PetscValidPointer(smooth, 3); 9889566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 989ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 99048a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 991ab8d36c9SPeter Brune *smooth = fas->smoothd; 9923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 993ab8d36c9SPeter Brune } 994ab8d36c9SPeter Brune 995ab8d36c9SPeter Brune /*@ 996ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level. 997ab8d36c9SPeter Brune 998ab8d36c9SPeter Brune Input Parameters: 999f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1000ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1001ab8d36c9SPeter Brune 1002f6dfbefdSBarry Smith Output Parameter: 1003ab8d36c9SPeter Brune smooth - the smoother 1004ab8d36c9SPeter Brune 1005ab8d36c9SPeter Brune Level: advanced 1006ab8d36c9SPeter Brune 1007f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1008ab8d36c9SPeter Brune @*/ 1009d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) 1010d71ae5a4SJacob Faibussowitsch { 1011ab8d36c9SPeter Brune SNES_FAS *fas; 1012ab8d36c9SPeter Brune SNES levelsnes; 101322d28d08SBarry Smith 1014ab8d36c9SPeter Brune PetscFunctionBegin; 1015f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1016f833ba53SLisandro Dalcin PetscValidPointer(smooth, 3); 10179566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 1018ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1019ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 102048a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 102148a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu)); 1022ab8d36c9SPeter Brune *smooth = fas->smoothd; 10233ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1024ab8d36c9SPeter Brune } 1025ab8d36c9SPeter Brune 1026ab8d36c9SPeter Brune /*@ 1027ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level. 1028ab8d36c9SPeter Brune 1029ab8d36c9SPeter Brune Input Parameters: 1030f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1031ab8d36c9SPeter Brune - level - the level (0 is coarsest) 1032ab8d36c9SPeter Brune 1033f6dfbefdSBarry Smith Output Parameter: 1034ab8d36c9SPeter Brune smooth - the smoother 1035ab8d36c9SPeter Brune 1036ab8d36c9SPeter Brune Level: advanced 1037ab8d36c9SPeter Brune 1038f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1039ab8d36c9SPeter Brune @*/ 1040d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) 1041d71ae5a4SJacob Faibussowitsch { 1042ab8d36c9SPeter Brune SNES_FAS *fas; 1043ab8d36c9SPeter Brune SNES levelsnes; 104422d28d08SBarry Smith 1045ab8d36c9SPeter Brune PetscFunctionBegin; 1046f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1047f833ba53SLisandro Dalcin PetscValidPointer(smooth, 3); 10489566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes)); 1049ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1050ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 105148a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 105248a46eb9SPierre Jolivet if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu)); 1053ab8d36c9SPeter Brune *smooth = fas->smoothu; 10543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1055ab8d36c9SPeter Brune } 1056d6ad1212SPeter Brune 1057d6ad1212SPeter Brune /*@ 1058d6ad1212SPeter Brune SNESFASGetCoarseSolve - Gets the coarsest solver. 1059d6ad1212SPeter Brune 1060f6dfbefdSBarry Smith Input Parameter: 1061f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 1062d6ad1212SPeter Brune 1063f6dfbefdSBarry Smith Output Parameter: 1064a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver 1065d6ad1212SPeter Brune 1066d6ad1212SPeter Brune Level: advanced 1067d6ad1212SPeter Brune 1068f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetInjection()`, `SNESFASSetRestriction()` 1069d6ad1212SPeter Brune @*/ 1070d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse) 1071d71ae5a4SJacob Faibussowitsch { 1072d6ad1212SPeter Brune SNES_FAS *fas; 1073d6ad1212SPeter Brune SNES levelsnes; 107422d28d08SBarry Smith 1075d6ad1212SPeter Brune PetscFunctionBegin; 1076f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1077064a246eSJacob Faibussowitsch PetscValidPointer(coarse, 2); 10789566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, 0, &levelsnes)); 1079d6ad1212SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1080d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 108148a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd)); 1082a3a80b83SMatthew G. Knepley *coarse = fas->smoothd; 10833ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1084d6ad1212SPeter Brune } 1085928e959bSPeter Brune 1086928e959bSPeter Brune /*@ 1087f6dfbefdSBarry Smith SNESFASFullSetDownSweep - Smooth during the initial downsweep for `SNESFAS` 1088928e959bSPeter Brune 1089c3339decSBarry Smith Logically Collective 1090928e959bSPeter Brune 1091928e959bSPeter Brune Input Parameters: 1092f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 1093928e959bSPeter Brune - swp - whether to downsweep or not 1094928e959bSPeter Brune 1095928e959bSPeter Brune Options Database Key: 1096928e959bSPeter Brune . -snes_fas_full_downsweep - Sets number of pre-smoothing steps 1097928e959bSPeter Brune 1098928e959bSPeter Brune Level: advanced 1099928e959bSPeter Brune 1100f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()` 1101928e959bSPeter Brune @*/ 1102d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullSetDownSweep(SNES snes, PetscBool swp) 1103d71ae5a4SJacob Faibussowitsch { 1104f833ba53SLisandro Dalcin SNES_FAS *fas; 1105928e959bSPeter Brune 1106928e959bSPeter Brune PetscFunctionBegin; 1107f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 1108f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 1109928e959bSPeter Brune fas->full_downsweep = swp; 11101baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetDownSweep(fas->next, swp)); 11113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 1112928e959bSPeter Brune } 11136dfbafefSToby Isaac 11146dfbafefSToby Isaac /*@ 11156dfbafefSToby Isaac SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 11166dfbafefSToby Isaac 1117c3339decSBarry Smith Logically Collective 11186dfbafefSToby Isaac 11196dfbafefSToby Isaac Input Parameters: 1120f6dfbefdSBarry Smith + snes - the `SNESFAS` nonlinear multigrid context 11216dfbafefSToby Isaac - total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 11226dfbafefSToby Isaac 11236dfbafefSToby Isaac Options Database Key: 112467b8a455SSatish Balay . -snes_fas_full_total - Use total restriction and interpolation on the initial down and up sweeps for the full FAS cycle 11256dfbafefSToby Isaac 11266dfbafefSToby Isaac Level: advanced 11276dfbafefSToby Isaac 1128f6dfbefdSBarry Smith Note: 1129f6dfbefdSBarry Smith This option is only significant if the interpolation of a coarse correction (`MatInterpolate()`) is significantly different from total 1130f6dfbefdSBarry Smith solution interpolation (`DMInterpolateSolution()`). 11316dfbafefSToby Isaac 1132f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()` 11336dfbafefSToby Isaac @*/ 1134d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullSetTotal(SNES snes, PetscBool total) 1135d71ae5a4SJacob Faibussowitsch { 11366dfbafefSToby Isaac SNES_FAS *fas; 11376dfbafefSToby Isaac 11386dfbafefSToby Isaac PetscFunctionBegin; 11396dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 11406dfbafefSToby Isaac fas = (SNES_FAS *)snes->data; 11416dfbafefSToby Isaac fas->full_total = total; 11421baa6e33SBarry Smith if (fas->next) PetscCall(SNESFASFullSetTotal(fas->next, total)); 11433ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11446dfbafefSToby Isaac } 11456dfbafefSToby Isaac 11466dfbafefSToby Isaac /*@ 11476dfbafefSToby Isaac SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 11486dfbafefSToby Isaac 1149c3339decSBarry Smith Logically Collective 11506dfbafefSToby Isaac 1151f6dfbefdSBarry Smith Input Parameter: 1152f6dfbefdSBarry Smith . snes - the `SNESFAS` nonlinear multigrid context 11536dfbafefSToby Isaac 11546dfbafefSToby Isaac Output: 11556dfbafefSToby Isaac . total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 11566dfbafefSToby Isaac 11576dfbafefSToby Isaac Level: advanced 11586dfbafefSToby Isaac 1159f6dfbefdSBarry Smith .seealso: `SNESFAS`, `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`, `SNESFullSetTotal()` 11606dfbafefSToby Isaac @*/ 1161d71ae5a4SJacob Faibussowitsch PetscErrorCode SNESFASFullGetTotal(SNES snes, PetscBool *total) 1162d71ae5a4SJacob Faibussowitsch { 11636dfbafefSToby Isaac SNES_FAS *fas; 11646dfbafefSToby Isaac 11656dfbafefSToby Isaac PetscFunctionBegin; 11666dfbafefSToby Isaac PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 11676dfbafefSToby Isaac fas = (SNES_FAS *)snes->data; 11686dfbafefSToby Isaac *total = fas->full_total; 11693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 11706dfbafefSToby Isaac } 1171