1a7b5fb5fSBarry Smith #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnes.h" I*/ 2ab8d36c9SPeter Brune 3ab8d36c9SPeter Brune /* -------------- functions called on the fine level -------------- */ 4ab8d36c9SPeter Brune 5ab8d36c9SPeter Brune /*@ 6ab8d36c9SPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 7ab8d36c9SPeter Brune 8ab8d36c9SPeter Brune Logically Collective 9ab8d36c9SPeter Brune 10ab8d36c9SPeter Brune Input Parameters: 11583a1250SSatish Balay + snes - FAS context 1234d65b3cSPeter Brune - fastype - SNES_FAS_ADDITIVE, SNES_FAS_MULTIPLICATIVE, SNES_FAS_FULL, or SNES_FAS_KASKADE 13583a1250SSatish Balay 14583a1250SSatish Balay Level: intermediate 15ab8d36c9SPeter Brune 16ab8d36c9SPeter Brune .seealso: PCMGSetType() 17ab8d36c9SPeter Brune @*/ 18ab8d36c9SPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 19ab8d36c9SPeter Brune { 20f833ba53SLisandro Dalcin SNES_FAS *fas; 21ab8d36c9SPeter Brune PetscErrorCode ierr; 2222d28d08SBarry Smith 23ab8d36c9SPeter Brune PetscFunctionBegin; 24f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 25ab8d36c9SPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 26f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 27ab8d36c9SPeter Brune fas->fastype = fastype; 2822d28d08SBarry Smith if (fas->next) { 2922d28d08SBarry Smith ierr = SNESFASSetType(fas->next, fastype);CHKERRQ(ierr); 3022d28d08SBarry Smith } 31ab8d36c9SPeter Brune PetscFunctionReturn(0); 32ab8d36c9SPeter Brune } 33ab8d36c9SPeter Brune 34ab8d36c9SPeter Brune /*@ 35ab8d36c9SPeter Brune SNESFASGetType - Sets the update and correction type used for FAS. 36ab8d36c9SPeter Brune 37ab8d36c9SPeter Brune Logically Collective 38ab8d36c9SPeter Brune 39ab8d36c9SPeter Brune Input Parameters: 40ab8d36c9SPeter Brune . snes - FAS context 41ab8d36c9SPeter Brune 42ab8d36c9SPeter Brune Output Parameters: 43ab8d36c9SPeter Brune . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE 44ab8d36c9SPeter Brune 45583a1250SSatish Balay Level: intermediate 46583a1250SSatish Balay 47ab8d36c9SPeter Brune .seealso: PCMGSetType() 48ab8d36c9SPeter Brune @*/ 49ab8d36c9SPeter Brune PetscErrorCode SNESFASGetType(SNES snes,SNESFASType *fastype) 50ab8d36c9SPeter Brune { 51f833ba53SLisandro Dalcin SNES_FAS *fas; 52ab8d36c9SPeter Brune 53ab8d36c9SPeter Brune PetscFunctionBegin; 54f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 55ab8d36c9SPeter Brune PetscValidPointer(fastype, 2); 56f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 57ab8d36c9SPeter Brune *fastype = fas->fastype; 58ab8d36c9SPeter Brune PetscFunctionReturn(0); 59ab8d36c9SPeter Brune } 60ab8d36c9SPeter Brune 61ab8d36c9SPeter Brune /*@C 62ab8d36c9SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 63ab8d36c9SPeter Brune Must be called before any other FAS routine. 64ab8d36c9SPeter Brune 65ab8d36c9SPeter Brune Input Parameters: 66ab8d36c9SPeter Brune + snes - the snes context 67ab8d36c9SPeter Brune . levels - the number of levels 68ab8d36c9SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 692bf68e3eSBarry Smith problems on smaller sets of processors. 70ab8d36c9SPeter Brune 71ab8d36c9SPeter Brune Level: intermediate 72ab8d36c9SPeter Brune 73ab8d36c9SPeter Brune Notes: 74ab8d36c9SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 75ab8d36c9SPeter Brune for setting the level options rather than the -fas_coarse prefix. 76ab8d36c9SPeter Brune 77ab8d36c9SPeter Brune .seealso: SNESFASGetLevels() 78ab8d36c9SPeter Brune @*/ 7922d28d08SBarry Smith PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms) 8022d28d08SBarry Smith { 81ab8d36c9SPeter Brune PetscErrorCode ierr; 82ab8d36c9SPeter Brune PetscInt i; 83ab8d36c9SPeter Brune const char *optionsprefix; 84ab8d36c9SPeter Brune char tprefix[128]; 85f833ba53SLisandro Dalcin SNES_FAS *fas; 86ab8d36c9SPeter Brune SNES prevsnes; 87ab8d36c9SPeter Brune MPI_Comm comm; 8822d28d08SBarry Smith 89ab8d36c9SPeter Brune PetscFunctionBegin; 90f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 91f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 92ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 93ab8d36c9SPeter Brune if (levels == fas->levels) { 9422d28d08SBarry Smith if (!comms) PetscFunctionReturn(0); 95ab8d36c9SPeter Brune } 96ab8d36c9SPeter Brune /* user has changed the number of levels; reset */ 97f833ba53SLisandro Dalcin ierr = (*snes->ops->reset)(snes);CHKERRQ(ierr); 98ab8d36c9SPeter Brune /* destroy any coarser levels if necessary */ 99f833ba53SLisandro Dalcin ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 1000298fd71SBarry Smith fas->next = NULL; 1010298fd71SBarry Smith fas->previous = NULL; 102ab8d36c9SPeter Brune prevsnes = snes; 103ab8d36c9SPeter Brune /* setup the finest level */ 104ab8d36c9SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 105ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1);CHKERRQ(ierr); 106ab8d36c9SPeter Brune for (i = levels - 1; i >= 0; i--) { 107ab8d36c9SPeter Brune if (comms) comm = comms[i]; 108ab8d36c9SPeter Brune fas->level = i; 109ab8d36c9SPeter Brune fas->levels = levels; 110ab8d36c9SPeter Brune fas->fine = snes; 1110298fd71SBarry Smith fas->next = NULL; 112ab8d36c9SPeter Brune if (i > 0) { 113ab8d36c9SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 114e964f0dbSPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 115f833ba53SLisandro Dalcin ierr = PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_cycle_",(int)fas->level);CHKERRQ(ierr); 116ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(fas->next,optionsprefix);CHKERRQ(ierr); 117e964f0dbSPeter Brune ierr = SNESAppendOptionsPrefix(fas->next,tprefix);CHKERRQ(ierr); 118ab8d36c9SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 119ab8d36c9SPeter Brune ierr = SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);CHKERRQ(ierr); 120ab8d36c9SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 121ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1);CHKERRQ(ierr); 1221aa26658SKarl Rupp 123ab8d36c9SPeter Brune ((SNES_FAS*)fas->next->data)->previous = prevsnes; 1241aa26658SKarl Rupp 125ab8d36c9SPeter Brune prevsnes = fas->next; 126ab8d36c9SPeter Brune fas = (SNES_FAS*)prevsnes->data; 127ab8d36c9SPeter Brune } 128ab8d36c9SPeter Brune } 129ab8d36c9SPeter Brune PetscFunctionReturn(0); 130ab8d36c9SPeter Brune } 131ab8d36c9SPeter Brune 132ab8d36c9SPeter Brune /*@ 133ab8d36c9SPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 134ab8d36c9SPeter Brune 135ab8d36c9SPeter Brune Input Parameter: 136ab8d36c9SPeter Brune . snes - the nonlinear solver context 137ab8d36c9SPeter Brune 138ab8d36c9SPeter Brune Output parameter: 139ab8d36c9SPeter Brune . levels - the number of levels 140ab8d36c9SPeter Brune 141ab8d36c9SPeter Brune Level: advanced 142ab8d36c9SPeter Brune 143ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 144ab8d36c9SPeter Brune @*/ 14522d28d08SBarry Smith PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels) 14622d28d08SBarry Smith { 147f833ba53SLisandro Dalcin SNES_FAS *fas; 1485fd66863SKarl Rupp 149ab8d36c9SPeter Brune PetscFunctionBegin; 150f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 151064a246eSJacob Faibussowitsch PetscValidIntPointer(levels,2); 152f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 153ab8d36c9SPeter Brune *levels = fas->levels; 154ab8d36c9SPeter Brune PetscFunctionReturn(0); 155ab8d36c9SPeter Brune } 156ab8d36c9SPeter Brune 157ab8d36c9SPeter Brune /*@ 158ab8d36c9SPeter Brune SNESFASGetCycleSNES - Gets the SNES corresponding to a particular 159ab8d36c9SPeter Brune level of the FAS hierarchy. 160ab8d36c9SPeter Brune 161ab8d36c9SPeter Brune Input Parameters: 162ab8d36c9SPeter Brune + snes - the multigrid context 163ab8d36c9SPeter Brune level - the level to get 164ab8d36c9SPeter Brune - lsnes - whether to use the nonlinear smoother or not 165ab8d36c9SPeter Brune 166ab8d36c9SPeter Brune Level: advanced 167ab8d36c9SPeter Brune 168ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 169ab8d36c9SPeter Brune @*/ 17022d28d08SBarry Smith PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes) 17122d28d08SBarry Smith { 172f833ba53SLisandro Dalcin SNES_FAS *fas; 173ab8d36c9SPeter Brune PetscInt i; 174ab8d36c9SPeter Brune 175ab8d36c9SPeter Brune PetscFunctionBegin; 176f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 177f833ba53SLisandro Dalcin PetscValidPointer(lsnes,3); 178f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 179*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(level > fas->levels-1,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels); 180*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(fas->level != fas->levels - 1,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level); 181ab8d36c9SPeter Brune 182ab8d36c9SPeter Brune *lsnes = snes; 183ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) { 184ab8d36c9SPeter Brune *lsnes = fas->next; 185ab8d36c9SPeter Brune fas = (SNES_FAS*)(*lsnes)->data; 186ab8d36c9SPeter Brune } 187*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(fas->level != level,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 188ab8d36c9SPeter Brune PetscFunctionReturn(0); 189ab8d36c9SPeter Brune } 190ab8d36c9SPeter Brune 191ab8d36c9SPeter Brune /*@ 192ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 193ab8d36c9SPeter Brune use on all levels. 194ab8d36c9SPeter Brune 195ab8d36c9SPeter Brune Logically Collective on SNES 196ab8d36c9SPeter Brune 197ab8d36c9SPeter Brune Input Parameters: 198ab8d36c9SPeter Brune + snes - the multigrid context 199ab8d36c9SPeter Brune - n - the number of smoothing steps 200ab8d36c9SPeter Brune 201ab8d36c9SPeter Brune Options Database Key: 202ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 203ab8d36c9SPeter Brune 204ab8d36c9SPeter Brune Level: advanced 205ab8d36c9SPeter Brune 206ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 207ab8d36c9SPeter Brune @*/ 20822d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) 20922d28d08SBarry Smith { 210f833ba53SLisandro Dalcin SNES_FAS *fas; 21122d28d08SBarry Smith PetscErrorCode ierr; 21222d28d08SBarry Smith 213ab8d36c9SPeter Brune PetscFunctionBegin; 214f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 215f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 216ab8d36c9SPeter Brune fas->max_up_it = n; 217656ede7eSPeter Brune if (!fas->smoothu && fas->level != 0) { 21822d28d08SBarry Smith ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 219ab8d36c9SPeter Brune } 22022d28d08SBarry Smith if (fas->smoothu) { 22122d28d08SBarry Smith ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr); 22222d28d08SBarry Smith } 223ab8d36c9SPeter Brune if (fas->next) { 224ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 225ab8d36c9SPeter Brune } 226ab8d36c9SPeter Brune PetscFunctionReturn(0); 227ab8d36c9SPeter Brune } 228ab8d36c9SPeter Brune 229ab8d36c9SPeter Brune /*@ 230ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 231ab8d36c9SPeter Brune use on all levels. 232ab8d36c9SPeter Brune 233ab8d36c9SPeter Brune Logically Collective on SNES 234ab8d36c9SPeter Brune 235ab8d36c9SPeter Brune Input Parameters: 236ab8d36c9SPeter Brune + snes - the multigrid context 237ab8d36c9SPeter Brune - n - the number of smoothing steps 238ab8d36c9SPeter Brune 239ab8d36c9SPeter Brune Options Database Key: 240ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 241ab8d36c9SPeter Brune 242ab8d36c9SPeter Brune Level: advanced 243ab8d36c9SPeter Brune 244ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 245ab8d36c9SPeter Brune @*/ 24622d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) 24722d28d08SBarry Smith { 248f833ba53SLisandro Dalcin SNES_FAS *fas; 249f833ba53SLisandro Dalcin PetscErrorCode ierr; 25022d28d08SBarry Smith 251ab8d36c9SPeter Brune PetscFunctionBegin; 252f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 253f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 254ab8d36c9SPeter Brune if (!fas->smoothd) { 25522d28d08SBarry Smith ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 256ab8d36c9SPeter Brune } 257ab8d36c9SPeter Brune ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr); 2581aa26658SKarl Rupp 259ab8d36c9SPeter Brune fas->max_down_it = n; 260ab8d36c9SPeter Brune if (fas->next) { 261ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 262ab8d36c9SPeter Brune } 263ab8d36c9SPeter Brune PetscFunctionReturn(0); 264ab8d36c9SPeter Brune } 265ab8d36c9SPeter Brune 26687f44e3fSPeter Brune /*@ 26787f44e3fSPeter Brune SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep 26887f44e3fSPeter Brune 26987f44e3fSPeter Brune Logically Collective on SNES 27087f44e3fSPeter Brune 27187f44e3fSPeter Brune Input Parameters: 27287f44e3fSPeter Brune + snes - the multigrid context 27387f44e3fSPeter Brune - n - the number of smoothing steps 27487f44e3fSPeter Brune 27587f44e3fSPeter Brune Options Database Key: 27687f44e3fSPeter Brune . -snes_fas_continuation - sets continuation to true 27787f44e3fSPeter Brune 27887f44e3fSPeter Brune Level: advanced 27987f44e3fSPeter Brune 28095452b02SPatrick Sanan Notes: 28195452b02SPatrick Sanan This sets the prefix on the upsweep smoothers to -fas_continuation 28287f44e3fSPeter Brune 28387f44e3fSPeter Brune .seealso: SNESFAS 28487f44e3fSPeter Brune @*/ 28587f44e3fSPeter Brune PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation) 28687f44e3fSPeter Brune { 28787f44e3fSPeter Brune const char *optionsprefix; 28887f44e3fSPeter Brune char tprefix[128]; 289f833ba53SLisandro Dalcin SNES_FAS *fas; 290f833ba53SLisandro Dalcin PetscErrorCode ierr; 29187f44e3fSPeter Brune 29287f44e3fSPeter Brune PetscFunctionBegin; 293f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 294f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 29587f44e3fSPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 29687f44e3fSPeter Brune if (!fas->smoothu) { 29787f44e3fSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 29887f44e3fSPeter Brune } 299f833ba53SLisandro Dalcin ierr = PetscStrncpy(tprefix,"fas_levels_continuation_",sizeof(tprefix));CHKERRQ(ierr); 30087f44e3fSPeter Brune ierr = SNESSetOptionsPrefix(fas->smoothu, optionsprefix);CHKERRQ(ierr); 30187f44e3fSPeter Brune ierr = SNESAppendOptionsPrefix(fas->smoothu, tprefix);CHKERRQ(ierr); 30287f44e3fSPeter Brune ierr = SNESSetType(fas->smoothu,SNESNEWTONLS);CHKERRQ(ierr); 30387f44e3fSPeter Brune ierr = SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100);CHKERRQ(ierr); 30487f44e3fSPeter Brune fas->continuation = continuation; 30587f44e3fSPeter Brune if (fas->next) { 30687f44e3fSPeter Brune ierr = SNESFASSetContinuation(fas->next,continuation);CHKERRQ(ierr); 30787f44e3fSPeter Brune } 30887f44e3fSPeter Brune PetscFunctionReturn(0); 30987f44e3fSPeter Brune } 31087f44e3fSPeter Brune 311ab8d36c9SPeter Brune /*@ 312ab8d36c9SPeter Brune SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 313ab8d36c9SPeter Brune complicated cycling. 314ab8d36c9SPeter Brune 315ab8d36c9SPeter Brune Logically Collective on SNES 316ab8d36c9SPeter Brune 317ab8d36c9SPeter Brune Input Parameters: 318ab8d36c9SPeter Brune + snes - the multigrid context 319ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 320ab8d36c9SPeter Brune 321ab8d36c9SPeter Brune Options Database Key: 322e1bc860dSBarry Smith . -snes_fas_cycles 1 or 2 323ab8d36c9SPeter Brune 324ab8d36c9SPeter Brune Level: advanced 325ab8d36c9SPeter Brune 326ab8d36c9SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 327ab8d36c9SPeter Brune @*/ 32822d28d08SBarry Smith PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) 32922d28d08SBarry Smith { 330f833ba53SLisandro Dalcin SNES_FAS *fas; 331ab8d36c9SPeter Brune PetscErrorCode ierr; 332ab8d36c9SPeter Brune PetscBool isFine; 33322d28d08SBarry Smith 334ab8d36c9SPeter Brune PetscFunctionBegin; 335f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 33622d28d08SBarry Smith ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 337f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 338ab8d36c9SPeter Brune fas->n_cycles = cycles; 3391aa26658SKarl Rupp if (!isFine) { 340ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 3411aa26658SKarl Rupp } 342ab8d36c9SPeter Brune if (fas->next) { 343ab8d36c9SPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 344ab8d36c9SPeter Brune } 345ab8d36c9SPeter Brune PetscFunctionReturn(0); 346ab8d36c9SPeter Brune } 347ab8d36c9SPeter Brune 348c8c899caSPeter Brune /*@ 349c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring 350c8c899caSPeter Brune 351c8c899caSPeter Brune Logically Collective on SNES 352c8c899caSPeter Brune 353c8c899caSPeter Brune Input Parameters: 354c8c899caSPeter Brune + snes - the FAS context 355d142ab34SLawrence Mitchell . vf - viewer and format structure (may be NULL if flg is FALSE) 356c8c899caSPeter Brune - flg - monitor or not 357c8c899caSPeter Brune 358c8c899caSPeter Brune Level: advanced 359c8c899caSPeter Brune 360c8c899caSPeter Brune .seealso: SNESFASSetCyclesOnLevel() 361c8c899caSPeter Brune @*/ 362d142ab34SLawrence Mitchell PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg) 36322d28d08SBarry Smith { 364f833ba53SLisandro Dalcin SNES_FAS *fas; 365c8c899caSPeter Brune PetscErrorCode ierr; 366c8c899caSPeter Brune PetscBool isFine; 367f833ba53SLisandro Dalcin PetscInt i, levels; 368c8c899caSPeter Brune SNES levelsnes; 36922d28d08SBarry Smith 370c8c899caSPeter Brune PetscFunctionBegin; 371f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 372c8c899caSPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 373f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 374f833ba53SLisandro Dalcin levels = fas->levels; 375c8c899caSPeter Brune if (isFine) { 376c8c899caSPeter Brune for (i = 0; i < levels; i++) { 37722d28d08SBarry Smith ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 378c8c899caSPeter Brune fas = (SNES_FAS*)levelsnes->data; 379c8c899caSPeter Brune if (flg) { 380c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */ 381c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 382d142ab34SLawrence Mitchell /* Only register destroy on finest level */ 383d142ab34SLawrence Mitchell ierr = SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL));CHKERRQ(ierr); 3841aa26658SKarl Rupp } else if (i != fas->levels - 1) { 385c8c899caSPeter Brune /* unset the monitors on the coarse levels */ 386c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 387c8c899caSPeter Brune } 388c8c899caSPeter Brune } 389c8c899caSPeter Brune } 390c8c899caSPeter Brune PetscFunctionReturn(0); 391c8c899caSPeter Brune } 392c8c899caSPeter Brune 3930dd27c6cSPeter Brune /*@ 3940dd27c6cSPeter Brune SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels 3950dd27c6cSPeter Brune 3960dd27c6cSPeter Brune Logically Collective on SNES 3970dd27c6cSPeter Brune 3980dd27c6cSPeter Brune Input Parameters: 3990dd27c6cSPeter Brune + snes - the FAS context 4000dd27c6cSPeter Brune - flg - monitor or not 4010dd27c6cSPeter Brune 4020dd27c6cSPeter Brune Level: advanced 4030dd27c6cSPeter Brune 4040dd27c6cSPeter Brune .seealso: SNESFASSetMonitor() 4050dd27c6cSPeter Brune @*/ 4060dd27c6cSPeter Brune PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) 4070dd27c6cSPeter Brune { 408f833ba53SLisandro Dalcin SNES_FAS *fas; 4090dd27c6cSPeter Brune PetscErrorCode ierr; 4100dd27c6cSPeter Brune PetscBool isFine; 411f833ba53SLisandro Dalcin PetscInt i, levels; 4120dd27c6cSPeter Brune SNES levelsnes; 4130dd27c6cSPeter Brune char eventname[128]; 4140dd27c6cSPeter Brune 4150dd27c6cSPeter Brune PetscFunctionBegin; 416f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 4170dd27c6cSPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 418f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 419f833ba53SLisandro Dalcin levels = fas->levels; 4200dd27c6cSPeter Brune if (isFine) { 4210dd27c6cSPeter Brune for (i = 0; i < levels; i++) { 4220dd27c6cSPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 4230dd27c6cSPeter Brune fas = (SNES_FAS*)levelsnes->data; 4240dd27c6cSPeter Brune if (flg) { 425f833ba53SLisandro Dalcin ierr = PetscSNPrintf(eventname,sizeof(eventname),"FASSetup %d",(int)i);CHKERRQ(ierr); 4260dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr); 427f833ba53SLisandro Dalcin ierr = PetscSNPrintf(eventname,sizeof(eventname),"FASSmooth %d",(int)i);CHKERRQ(ierr); 4280dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr); 429f833ba53SLisandro Dalcin ierr = PetscSNPrintf(eventname,sizeof(eventname),"FASResid %d",(int)i);CHKERRQ(ierr); 4300dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr); 431f833ba53SLisandro Dalcin ierr = PetscSNPrintf(eventname,sizeof(eventname),"FASInterp %d",(int)i);CHKERRQ(ierr); 4320dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr); 4330dd27c6cSPeter Brune } else { 4340298fd71SBarry Smith fas->eventsmoothsetup = 0; 4350298fd71SBarry Smith fas->eventsmoothsolve = 0; 4360298fd71SBarry Smith fas->eventresidual = 0; 4370298fd71SBarry Smith fas->eventinterprestrict = 0; 4380dd27c6cSPeter Brune } 4390dd27c6cSPeter Brune } 4400dd27c6cSPeter Brune } 4410dd27c6cSPeter Brune PetscFunctionReturn(0); 4420dd27c6cSPeter Brune } 4430dd27c6cSPeter Brune 444ab8d36c9SPeter Brune /* 445ab8d36c9SPeter Brune Creates the default smoother type. 446ab8d36c9SPeter Brune 44704d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level. 448ab8d36c9SPeter Brune 449ab8d36c9SPeter Brune */ 45022d28d08SBarry Smith PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) 45122d28d08SBarry Smith { 452ab8d36c9SPeter Brune SNES_FAS *fas; 453ab8d36c9SPeter Brune const char *optionsprefix; 454ab8d36c9SPeter Brune char tprefix[128]; 455ab8d36c9SPeter Brune PetscErrorCode ierr; 456ab8d36c9SPeter Brune SNES nsmooth; 45722d28d08SBarry Smith 458ab8d36c9SPeter Brune PetscFunctionBegin; 459f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 460064a246eSJacob Faibussowitsch PetscValidPointer(smooth,2); 461ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 462ab8d36c9SPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 463ab8d36c9SPeter Brune /* create the default smoother */ 464ce94432eSBarry Smith ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);CHKERRQ(ierr); 465ab8d36c9SPeter Brune if (fas->level == 0) { 466f833ba53SLisandro Dalcin ierr = PetscStrncpy(tprefix,"fas_coarse_",sizeof(tprefix));CHKERRQ(ierr); 467ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 468ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 46904d7464bSBarry Smith ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr); 470e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr); 471ab8d36c9SPeter Brune } else { 472f833ba53SLisandro Dalcin ierr = PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_",(int)fas->level);CHKERRQ(ierr); 473ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 474ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 475ab8d36c9SPeter Brune ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr); 476e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr); 477ab8d36c9SPeter Brune } 478ab8d36c9SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 4793bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);CHKERRQ(ierr); 480f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr); 481ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level);CHKERRQ(ierr); 482ab8d36c9SPeter Brune *smooth = nsmooth; 483ab8d36c9SPeter Brune PetscFunctionReturn(0); 484ab8d36c9SPeter Brune } 485ab8d36c9SPeter Brune 486ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */ 487ab8d36c9SPeter Brune 488ab8d36c9SPeter Brune /*@ 489ab8d36c9SPeter Brune SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 490ab8d36c9SPeter Brune 491ab8d36c9SPeter Brune Logically Collective on SNES 492ab8d36c9SPeter Brune 493ab8d36c9SPeter Brune Input Parameters: 494ab8d36c9SPeter Brune + snes - the multigrid context 495ab8d36c9SPeter Brune . level - the level to set the number of cycles on 496ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 497ab8d36c9SPeter Brune 498ab8d36c9SPeter Brune Level: advanced 499ab8d36c9SPeter Brune 500ab8d36c9SPeter Brune .seealso: SNESFASSetCycles() 501ab8d36c9SPeter Brune @*/ 50222d28d08SBarry Smith PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) 50322d28d08SBarry Smith { 504f833ba53SLisandro Dalcin SNES_FAS *fas; 505ab8d36c9SPeter Brune PetscErrorCode ierr; 50622d28d08SBarry Smith 507ab8d36c9SPeter Brune PetscFunctionBegin; 508f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 509f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 510ab8d36c9SPeter Brune fas->n_cycles = cycles; 511ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 512ab8d36c9SPeter Brune PetscFunctionReturn(0); 513ab8d36c9SPeter Brune } 514ab8d36c9SPeter Brune 515ab8d36c9SPeter Brune /*@ 516ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 517ab8d36c9SPeter Brune 518ab8d36c9SPeter Brune Logically Collective on SNES 519ab8d36c9SPeter Brune 520ab8d36c9SPeter Brune Input Parameters: 521ab8d36c9SPeter Brune . snes - the multigrid context 522ab8d36c9SPeter Brune 523ab8d36c9SPeter Brune Output Parameters: 524ab8d36c9SPeter Brune . smooth - the smoother 525ab8d36c9SPeter Brune 526ab8d36c9SPeter Brune Level: advanced 527ab8d36c9SPeter Brune 528ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown() 529ab8d36c9SPeter Brune @*/ 530ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 531ab8d36c9SPeter Brune { 532ab8d36c9SPeter Brune SNES_FAS *fas; 5335fd66863SKarl Rupp 534ab8d36c9SPeter Brune PetscFunctionBegin; 535f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 536f833ba53SLisandro Dalcin PetscValidPointer(smooth,2); 537ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 538ab8d36c9SPeter Brune *smooth = fas->smoothd; 539ab8d36c9SPeter Brune PetscFunctionReturn(0); 540ab8d36c9SPeter Brune } 541ab8d36c9SPeter Brune /*@ 542ab8d36c9SPeter Brune SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 543ab8d36c9SPeter Brune 544ab8d36c9SPeter Brune Logically Collective on SNES 545ab8d36c9SPeter Brune 546ab8d36c9SPeter Brune Input Parameters: 547ab8d36c9SPeter Brune . snes - the multigrid context 548ab8d36c9SPeter Brune 549ab8d36c9SPeter Brune Output Parameters: 550ab8d36c9SPeter Brune . smoothu - the smoother 551ab8d36c9SPeter Brune 552ab8d36c9SPeter Brune Notes: 553ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent 554ab8d36c9SPeter Brune default behavior in the process of the solve. 555ab8d36c9SPeter Brune 556ab8d36c9SPeter Brune Level: advanced 557ab8d36c9SPeter Brune 558ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown() 559ab8d36c9SPeter Brune @*/ 560ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 561ab8d36c9SPeter Brune { 562ab8d36c9SPeter Brune SNES_FAS *fas; 5635fd66863SKarl Rupp 564ab8d36c9SPeter Brune PetscFunctionBegin; 565f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 566f833ba53SLisandro Dalcin PetscValidPointer(smoothu,2); 567ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 5681aa26658SKarl Rupp if (!fas->smoothu) *smoothu = fas->smoothd; 5691aa26658SKarl Rupp else *smoothu = fas->smoothu; 570ab8d36c9SPeter Brune PetscFunctionReturn(0); 571ab8d36c9SPeter Brune } 572ab8d36c9SPeter Brune 573ab8d36c9SPeter Brune /*@ 574ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 575ab8d36c9SPeter Brune 576ab8d36c9SPeter Brune Logically Collective on SNES 577ab8d36c9SPeter Brune 578ab8d36c9SPeter Brune Input Parameters: 579ab8d36c9SPeter Brune . snes - the multigrid context 580ab8d36c9SPeter Brune 581ab8d36c9SPeter Brune Output Parameters: 582ab8d36c9SPeter Brune . smoothd - the smoother 583ab8d36c9SPeter Brune 584ab8d36c9SPeter Brune Level: advanced 585ab8d36c9SPeter Brune 586ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 587ab8d36c9SPeter Brune @*/ 588ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 589ab8d36c9SPeter Brune { 590ab8d36c9SPeter Brune SNES_FAS *fas; 5915fd66863SKarl Rupp 592ab8d36c9SPeter Brune PetscFunctionBegin; 593f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 594f833ba53SLisandro Dalcin PetscValidPointer(smoothd,2); 595ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 596ab8d36c9SPeter Brune *smoothd = fas->smoothd; 597ab8d36c9SPeter Brune PetscFunctionReturn(0); 598ab8d36c9SPeter Brune } 599ab8d36c9SPeter Brune 600ab8d36c9SPeter Brune /*@ 601ab8d36c9SPeter Brune SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 602ab8d36c9SPeter Brune 603ab8d36c9SPeter Brune Logically Collective on SNES 604ab8d36c9SPeter Brune 605ab8d36c9SPeter Brune Input Parameters: 606ab8d36c9SPeter Brune . snes - the multigrid context 607ab8d36c9SPeter Brune 608ab8d36c9SPeter Brune Output Parameters: 609ab8d36c9SPeter Brune . correction - the coarse correction on this level 610ab8d36c9SPeter Brune 611ab8d36c9SPeter Brune Notes: 6120298fd71SBarry Smith Returns NULL on the coarsest level. 613ab8d36c9SPeter Brune 614ab8d36c9SPeter Brune Level: advanced 615ab8d36c9SPeter Brune 616ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 617ab8d36c9SPeter Brune @*/ 618ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 619ab8d36c9SPeter Brune { 620ab8d36c9SPeter Brune SNES_FAS *fas; 6215fd66863SKarl Rupp 622ab8d36c9SPeter Brune PetscFunctionBegin; 623f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 624f833ba53SLisandro Dalcin PetscValidPointer(correction,2); 625ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 626ab8d36c9SPeter Brune *correction = fas->next; 627ab8d36c9SPeter Brune PetscFunctionReturn(0); 628ab8d36c9SPeter Brune } 629ab8d36c9SPeter Brune 630ab8d36c9SPeter Brune /*@ 631ab8d36c9SPeter Brune SNESFASCycleGetInterpolation - Gets the interpolation on this level 632ab8d36c9SPeter Brune 633ab8d36c9SPeter Brune Logically Collective on SNES 634ab8d36c9SPeter Brune 635ab8d36c9SPeter Brune Input Parameters: 636ab8d36c9SPeter Brune . snes - the multigrid context 637ab8d36c9SPeter Brune 638ab8d36c9SPeter Brune Output Parameters: 639ab8d36c9SPeter Brune . mat - the interpolation operator on this level 640ab8d36c9SPeter Brune 641ab8d36c9SPeter Brune Level: developer 642ab8d36c9SPeter Brune 643ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 644ab8d36c9SPeter Brune @*/ 645ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 646ab8d36c9SPeter Brune { 647ab8d36c9SPeter Brune SNES_FAS *fas; 6485fd66863SKarl Rupp 649ab8d36c9SPeter Brune PetscFunctionBegin; 650f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 651f833ba53SLisandro Dalcin PetscValidPointer(mat,2); 652ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 653ab8d36c9SPeter Brune *mat = fas->interpolate; 654ab8d36c9SPeter Brune PetscFunctionReturn(0); 655ab8d36c9SPeter Brune } 656ab8d36c9SPeter Brune 657ab8d36c9SPeter Brune /*@ 658ab8d36c9SPeter Brune SNESFASCycleGetRestriction - Gets the restriction on this level 659ab8d36c9SPeter Brune 660ab8d36c9SPeter Brune Logically Collective on SNES 661ab8d36c9SPeter Brune 662ab8d36c9SPeter Brune Input Parameters: 663ab8d36c9SPeter Brune . snes - the multigrid context 664ab8d36c9SPeter Brune 665ab8d36c9SPeter Brune Output Parameters: 666ab8d36c9SPeter Brune . mat - the restriction operator on this level 667ab8d36c9SPeter Brune 668ab8d36c9SPeter Brune Level: developer 669ab8d36c9SPeter Brune 670ab8d36c9SPeter Brune .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation() 671ab8d36c9SPeter Brune @*/ 672ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 673ab8d36c9SPeter Brune { 674ab8d36c9SPeter Brune SNES_FAS *fas; 6755fd66863SKarl Rupp 676ab8d36c9SPeter Brune PetscFunctionBegin; 677f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 678f833ba53SLisandro Dalcin PetscValidPointer(mat,2); 679ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 680ab8d36c9SPeter Brune *mat = fas->restrct; 681ab8d36c9SPeter Brune PetscFunctionReturn(0); 682ab8d36c9SPeter Brune } 683ab8d36c9SPeter Brune 684ab8d36c9SPeter Brune /*@ 685ab8d36c9SPeter Brune SNESFASCycleGetInjection - Gets the injection on this level 686ab8d36c9SPeter Brune 687ab8d36c9SPeter Brune Logically Collective on SNES 688ab8d36c9SPeter Brune 689ab8d36c9SPeter Brune Input Parameters: 690ab8d36c9SPeter Brune . snes - the multigrid context 691ab8d36c9SPeter Brune 692ab8d36c9SPeter Brune Output Parameters: 693ab8d36c9SPeter Brune . mat - the restriction operator on this level 694ab8d36c9SPeter Brune 695ab8d36c9SPeter Brune Level: developer 696ab8d36c9SPeter Brune 697ab8d36c9SPeter Brune .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction() 698ab8d36c9SPeter Brune @*/ 699ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 700ab8d36c9SPeter Brune { 701ab8d36c9SPeter Brune SNES_FAS *fas; 7025fd66863SKarl Rupp 703ab8d36c9SPeter Brune PetscFunctionBegin; 704f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 705f833ba53SLisandro Dalcin PetscValidPointer(mat,2); 706ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 707ab8d36c9SPeter Brune *mat = fas->inject; 708ab8d36c9SPeter Brune PetscFunctionReturn(0); 709ab8d36c9SPeter Brune } 710ab8d36c9SPeter Brune 711ab8d36c9SPeter Brune /*@ 712ab8d36c9SPeter Brune SNESFASCycleGetRScale - Gets the injection on this level 713ab8d36c9SPeter Brune 714ab8d36c9SPeter Brune Logically Collective on SNES 715ab8d36c9SPeter Brune 716ab8d36c9SPeter Brune Input Parameters: 717ab8d36c9SPeter Brune . snes - the multigrid context 718ab8d36c9SPeter Brune 719ab8d36c9SPeter Brune Output Parameters: 720ab8d36c9SPeter Brune . mat - the restriction operator on this level 721ab8d36c9SPeter Brune 722ab8d36c9SPeter Brune Level: developer 723ab8d36c9SPeter Brune 724ab8d36c9SPeter Brune .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale() 725ab8d36c9SPeter Brune @*/ 726ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 727ab8d36c9SPeter Brune { 728ab8d36c9SPeter Brune SNES_FAS *fas; 7295fd66863SKarl Rupp 730ab8d36c9SPeter Brune PetscFunctionBegin; 731f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 732f833ba53SLisandro Dalcin PetscValidPointer(vec,2); 733ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 734ab8d36c9SPeter Brune *vec = fas->rscale; 735ab8d36c9SPeter Brune PetscFunctionReturn(0); 736ab8d36c9SPeter Brune } 737ab8d36c9SPeter Brune 738ab8d36c9SPeter Brune /*@ 739ab8d36c9SPeter Brune SNESFASCycleIsFine - Determines if a given cycle is the fine level. 740ab8d36c9SPeter Brune 741ab8d36c9SPeter Brune Logically Collective on SNES 742ab8d36c9SPeter Brune 743ab8d36c9SPeter Brune Input Parameters: 744ab8d36c9SPeter Brune . snes - the FAS context 745ab8d36c9SPeter Brune 746ab8d36c9SPeter Brune Output Parameters: 747ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not 748ab8d36c9SPeter Brune 749ab8d36c9SPeter Brune Level: advanced 750ab8d36c9SPeter Brune 751ab8d36c9SPeter Brune .seealso: SNESFASSetLevels() 752ab8d36c9SPeter Brune @*/ 753ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 754ab8d36c9SPeter Brune { 755ab8d36c9SPeter Brune SNES_FAS *fas; 7565fd66863SKarl Rupp 757ab8d36c9SPeter Brune PetscFunctionBegin; 758f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 759534a8f05SLisandro Dalcin PetscValidBoolPointer(flg,2); 760ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 7611aa26658SKarl Rupp if (fas->level == fas->levels - 1) *flg = PETSC_TRUE; 7621aa26658SKarl Rupp else *flg = PETSC_FALSE; 763ab8d36c9SPeter Brune PetscFunctionReturn(0); 764ab8d36c9SPeter Brune } 765ab8d36c9SPeter Brune 766ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */ 767ab8d36c9SPeter Brune 768ab8d36c9SPeter Brune /*@ 769ab8d36c9SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 770ab8d36c9SPeter Brune interpolation from l-1 to the lth level 771ab8d36c9SPeter Brune 772ab8d36c9SPeter Brune Input Parameters: 773ab8d36c9SPeter Brune + snes - the multigrid context 774ab8d36c9SPeter Brune . mat - the interpolation operator 775ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 776ab8d36c9SPeter Brune 777ab8d36c9SPeter Brune Level: advanced 778ab8d36c9SPeter Brune 779ab8d36c9SPeter Brune Notes: 780ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction 781ab8d36c9SPeter Brune for the same level. 782ab8d36c9SPeter Brune 783ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 784ab8d36c9SPeter Brune out from the matrix size which one it is. 785ab8d36c9SPeter Brune 786ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 787ab8d36c9SPeter Brune @*/ 7880adebc6cSBarry Smith PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) 7890adebc6cSBarry Smith { 79022d28d08SBarry Smith SNES_FAS *fas; 791ab8d36c9SPeter Brune PetscErrorCode ierr; 792ab8d36c9SPeter Brune SNES levelsnes; 79322d28d08SBarry Smith 794ab8d36c9SPeter Brune PetscFunctionBegin; 795f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 796f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 797ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 798ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 799ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 800ab8d36c9SPeter Brune ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 801ab8d36c9SPeter Brune fas->interpolate = mat; 802ab8d36c9SPeter Brune PetscFunctionReturn(0); 803ab8d36c9SPeter Brune } 804ab8d36c9SPeter Brune 805ab8d36c9SPeter Brune /*@ 806ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the 807ab8d36c9SPeter Brune interpolation from l-1 to the lth level 808ab8d36c9SPeter Brune 809ab8d36c9SPeter Brune Input Parameters: 810ab8d36c9SPeter Brune + snes - the multigrid context 811ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 812ab8d36c9SPeter Brune 813ab8d36c9SPeter Brune Output Parameters: 814ab8d36c9SPeter Brune . mat - the interpolation operator 815ab8d36c9SPeter Brune 816ab8d36c9SPeter Brune Level: advanced 817ab8d36c9SPeter Brune 818ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale() 819ab8d36c9SPeter Brune @*/ 82022d28d08SBarry Smith PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) 82122d28d08SBarry Smith { 82222d28d08SBarry Smith SNES_FAS *fas; 823ab8d36c9SPeter Brune PetscErrorCode ierr; 824ab8d36c9SPeter Brune SNES levelsnes; 82522d28d08SBarry Smith 826ab8d36c9SPeter Brune PetscFunctionBegin; 827f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 828f833ba53SLisandro Dalcin PetscValidPointer(mat,3); 829ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 830ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 831ab8d36c9SPeter Brune *mat = fas->interpolate; 832ab8d36c9SPeter Brune PetscFunctionReturn(0); 833ab8d36c9SPeter Brune } 834ab8d36c9SPeter Brune 835ab8d36c9SPeter Brune /*@ 836ab8d36c9SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 837ab8d36c9SPeter Brune from level l to l-1. 838ab8d36c9SPeter Brune 839ab8d36c9SPeter Brune Input Parameters: 840ab8d36c9SPeter Brune + snes - the multigrid context 841ab8d36c9SPeter Brune . mat - the restriction matrix 842ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 843ab8d36c9SPeter Brune 844ab8d36c9SPeter Brune Level: advanced 845ab8d36c9SPeter Brune 846ab8d36c9SPeter Brune Notes: 847ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation 848ab8d36c9SPeter Brune for the same level. 849ab8d36c9SPeter Brune 850ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 851ab8d36c9SPeter Brune out from the matrix size which one it is. 852ab8d36c9SPeter Brune 853ab8d36c9SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 854ab8d36c9SPeter Brune is used. 855ab8d36c9SPeter Brune 856ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 857ab8d36c9SPeter Brune @*/ 85822d28d08SBarry Smith PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) 85922d28d08SBarry Smith { 86022d28d08SBarry Smith SNES_FAS *fas; 861ab8d36c9SPeter Brune PetscErrorCode ierr; 862ab8d36c9SPeter Brune SNES levelsnes; 86322d28d08SBarry Smith 864ab8d36c9SPeter Brune PetscFunctionBegin; 865f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 866f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 867ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 868ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 869ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 870ab8d36c9SPeter Brune ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 871ab8d36c9SPeter Brune fas->restrct = mat; 872ab8d36c9SPeter Brune PetscFunctionReturn(0); 873ab8d36c9SPeter Brune } 874ab8d36c9SPeter Brune 875ab8d36c9SPeter Brune /*@ 876ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the 877ab8d36c9SPeter Brune restriction from l to the l-1th level 878ab8d36c9SPeter Brune 879ab8d36c9SPeter Brune Input Parameters: 880ab8d36c9SPeter Brune + snes - the multigrid context 881ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 882ab8d36c9SPeter Brune 883ab8d36c9SPeter Brune Output Parameters: 884ab8d36c9SPeter Brune . mat - the interpolation operator 885ab8d36c9SPeter Brune 886ab8d36c9SPeter Brune Level: advanced 887ab8d36c9SPeter Brune 888ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale() 889ab8d36c9SPeter Brune @*/ 89022d28d08SBarry Smith PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) 89122d28d08SBarry Smith { 89222d28d08SBarry Smith SNES_FAS *fas; 893ab8d36c9SPeter Brune PetscErrorCode ierr; 894ab8d36c9SPeter Brune SNES levelsnes; 89522d28d08SBarry Smith 896ab8d36c9SPeter Brune PetscFunctionBegin; 897f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 898f833ba53SLisandro Dalcin PetscValidPointer(mat,3); 899ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 900ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 901ab8d36c9SPeter Brune *mat = fas->restrct; 902ab8d36c9SPeter Brune PetscFunctionReturn(0); 903ab8d36c9SPeter Brune } 904ab8d36c9SPeter Brune 905ab8d36c9SPeter Brune /*@ 906ab8d36c9SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 907ab8d36c9SPeter Brune from level l to l-1. 908ab8d36c9SPeter Brune 909ab8d36c9SPeter Brune Input Parameters: 910ab8d36c9SPeter Brune + snes - the multigrid context 911ab8d36c9SPeter Brune . mat - the restriction matrix 912ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 913ab8d36c9SPeter Brune 914ab8d36c9SPeter Brune Level: advanced 915ab8d36c9SPeter Brune 916ab8d36c9SPeter Brune Notes: 917ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to 918ab8d36c9SPeter Brune project the solution instead. 919ab8d36c9SPeter Brune 920ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 921ab8d36c9SPeter Brune @*/ 92222d28d08SBarry Smith PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) 92322d28d08SBarry Smith { 92422d28d08SBarry Smith SNES_FAS *fas; 925ab8d36c9SPeter Brune PetscErrorCode ierr; 926ab8d36c9SPeter Brune SNES levelsnes; 92722d28d08SBarry Smith 928ab8d36c9SPeter Brune PetscFunctionBegin; 929f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 930f833ba53SLisandro Dalcin if (mat) PetscValidHeaderSpecific(mat,MAT_CLASSID,3); 931ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 932ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 933ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 934ab8d36c9SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 9351aa26658SKarl Rupp 936ab8d36c9SPeter Brune fas->inject = mat; 937ab8d36c9SPeter Brune PetscFunctionReturn(0); 938ab8d36c9SPeter Brune } 939ab8d36c9SPeter Brune 940ab8d36c9SPeter Brune /*@ 941ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the 942ab8d36c9SPeter Brune injection from l-1 to the lth level 943ab8d36c9SPeter Brune 944ab8d36c9SPeter Brune Input Parameters: 945ab8d36c9SPeter Brune + snes - the multigrid context 946ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 947ab8d36c9SPeter Brune 948ab8d36c9SPeter Brune Output Parameters: 949ab8d36c9SPeter Brune . mat - the injection operator 950ab8d36c9SPeter Brune 951ab8d36c9SPeter Brune Level: advanced 952ab8d36c9SPeter Brune 953ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale() 954ab8d36c9SPeter Brune @*/ 95522d28d08SBarry Smith PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) 95622d28d08SBarry Smith { 95722d28d08SBarry Smith SNES_FAS *fas; 958ab8d36c9SPeter Brune PetscErrorCode ierr; 959ab8d36c9SPeter Brune SNES levelsnes; 96022d28d08SBarry Smith 961ab8d36c9SPeter Brune PetscFunctionBegin; 962f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 963f833ba53SLisandro Dalcin PetscValidPointer(mat,3); 964ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 965ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 966ab8d36c9SPeter Brune *mat = fas->inject; 967ab8d36c9SPeter Brune PetscFunctionReturn(0); 968ab8d36c9SPeter Brune } 969ab8d36c9SPeter Brune 970ab8d36c9SPeter Brune /*@ 971ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 972ab8d36c9SPeter Brune operator from level l to l-1. 973ab8d36c9SPeter Brune 974ab8d36c9SPeter Brune Input Parameters: 975ab8d36c9SPeter Brune + snes - the multigrid context 976ab8d36c9SPeter Brune . rscale - the restriction scaling 977ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 978ab8d36c9SPeter Brune 979ab8d36c9SPeter Brune Level: advanced 980ab8d36c9SPeter Brune 981ab8d36c9SPeter Brune Notes: 982ab8d36c9SPeter Brune This is only used in the case that the injection is not set. 983ab8d36c9SPeter Brune 984ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 985ab8d36c9SPeter Brune @*/ 98622d28d08SBarry Smith PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) 98722d28d08SBarry Smith { 988ab8d36c9SPeter Brune SNES_FAS *fas; 989ab8d36c9SPeter Brune PetscErrorCode ierr; 990ab8d36c9SPeter Brune SNES levelsnes; 99122d28d08SBarry Smith 992ab8d36c9SPeter Brune PetscFunctionBegin; 993f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 994f833ba53SLisandro Dalcin if (rscale) PetscValidHeaderSpecific(rscale,VEC_CLASSID,3); 995ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 996ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 997ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 998ab8d36c9SPeter Brune ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 999ab8d36c9SPeter Brune fas->rscale = rscale; 1000ab8d36c9SPeter Brune PetscFunctionReturn(0); 1001ab8d36c9SPeter Brune } 1002ab8d36c9SPeter Brune 1003ab8d36c9SPeter Brune /*@ 1004ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level. 1005ab8d36c9SPeter Brune 1006ab8d36c9SPeter Brune Input Parameters: 1007ab8d36c9SPeter Brune + snes - the multigrid context 1008ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1009ab8d36c9SPeter Brune 1010ab8d36c9SPeter Brune Output Parameters: 1011ab8d36c9SPeter Brune smooth - the smoother 1012ab8d36c9SPeter Brune 1013ab8d36c9SPeter Brune Level: advanced 1014ab8d36c9SPeter Brune 1015ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1016ab8d36c9SPeter Brune @*/ 101722d28d08SBarry Smith PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) 101822d28d08SBarry Smith { 1019ab8d36c9SPeter Brune SNES_FAS *fas; 1020ab8d36c9SPeter Brune PetscErrorCode ierr; 1021ab8d36c9SPeter Brune SNES levelsnes; 102222d28d08SBarry Smith 1023ab8d36c9SPeter Brune PetscFunctionBegin; 1024f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 1025f833ba53SLisandro Dalcin PetscValidPointer(smooth,3); 1026ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1027ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1028ab8d36c9SPeter Brune if (!fas->smoothd) { 10293de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1030ab8d36c9SPeter Brune } 1031ab8d36c9SPeter Brune *smooth = fas->smoothd; 1032ab8d36c9SPeter Brune PetscFunctionReturn(0); 1033ab8d36c9SPeter Brune } 1034ab8d36c9SPeter Brune 1035ab8d36c9SPeter Brune /*@ 1036ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level. 1037ab8d36c9SPeter Brune 1038ab8d36c9SPeter Brune Input Parameters: 1039ab8d36c9SPeter Brune + snes - the multigrid context 1040ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1041ab8d36c9SPeter Brune 1042ab8d36c9SPeter Brune Output Parameters: 1043ab8d36c9SPeter Brune smooth - the smoother 1044ab8d36c9SPeter Brune 1045ab8d36c9SPeter Brune Level: advanced 1046ab8d36c9SPeter Brune 1047ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1048ab8d36c9SPeter Brune @*/ 104922d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) 105022d28d08SBarry Smith { 1051ab8d36c9SPeter Brune SNES_FAS *fas; 1052ab8d36c9SPeter Brune PetscErrorCode ierr; 1053ab8d36c9SPeter Brune SNES levelsnes; 105422d28d08SBarry Smith 1055ab8d36c9SPeter Brune PetscFunctionBegin; 1056f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 1057f833ba53SLisandro Dalcin PetscValidPointer(smooth,3); 1058ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1059ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1060ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1061ab8d36c9SPeter Brune if (!fas->smoothd) { 10623de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1063ab8d36c9SPeter Brune } 1064ab8d36c9SPeter Brune if (!fas->smoothu) { 10653de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr); 1066ab8d36c9SPeter Brune } 1067ab8d36c9SPeter Brune *smooth = fas->smoothd; 1068ab8d36c9SPeter Brune PetscFunctionReturn(0); 1069ab8d36c9SPeter Brune } 1070ab8d36c9SPeter Brune 1071ab8d36c9SPeter Brune /*@ 1072ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level. 1073ab8d36c9SPeter Brune 1074ab8d36c9SPeter Brune Input Parameters: 1075ab8d36c9SPeter Brune + snes - the multigrid context 1076ab8d36c9SPeter Brune - level - the level (0 is coarsest) 1077ab8d36c9SPeter Brune 1078ab8d36c9SPeter Brune Output Parameters: 1079ab8d36c9SPeter Brune smooth - the smoother 1080ab8d36c9SPeter Brune 1081ab8d36c9SPeter Brune Level: advanced 1082ab8d36c9SPeter Brune 1083ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1084ab8d36c9SPeter Brune @*/ 108522d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) 108622d28d08SBarry Smith { 1087ab8d36c9SPeter Brune SNES_FAS *fas; 1088ab8d36c9SPeter Brune PetscErrorCode ierr; 1089ab8d36c9SPeter Brune SNES levelsnes; 109022d28d08SBarry Smith 1091ab8d36c9SPeter Brune PetscFunctionBegin; 1092f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 1093f833ba53SLisandro Dalcin PetscValidPointer(smooth,3); 1094ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1095ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1096ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1097ab8d36c9SPeter Brune if (!fas->smoothd) { 10983de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1099ab8d36c9SPeter Brune } 1100ab8d36c9SPeter Brune if (!fas->smoothu) { 11013de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr); 1102ab8d36c9SPeter Brune } 1103ab8d36c9SPeter Brune *smooth = fas->smoothu; 1104ab8d36c9SPeter Brune PetscFunctionReturn(0); 1105ab8d36c9SPeter Brune } 1106d6ad1212SPeter Brune 1107d6ad1212SPeter Brune /*@ 1108d6ad1212SPeter Brune SNESFASGetCoarseSolve - Gets the coarsest solver. 1109d6ad1212SPeter Brune 1110d6ad1212SPeter Brune Input Parameters: 1111a3a80b83SMatthew G. Knepley . snes - the multigrid context 1112d6ad1212SPeter Brune 1113d6ad1212SPeter Brune Output Parameters: 1114a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver 1115d6ad1212SPeter Brune 1116d6ad1212SPeter Brune Level: advanced 1117d6ad1212SPeter Brune 1118d6ad1212SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1119d6ad1212SPeter Brune @*/ 1120a3a80b83SMatthew G. Knepley PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse) 112122d28d08SBarry Smith { 1122d6ad1212SPeter Brune SNES_FAS *fas; 1123d6ad1212SPeter Brune PetscErrorCode ierr; 1124d6ad1212SPeter Brune SNES levelsnes; 112522d28d08SBarry Smith 1126d6ad1212SPeter Brune PetscFunctionBegin; 1127f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 1128064a246eSJacob Faibussowitsch PetscValidPointer(coarse,2); 1129d6ad1212SPeter Brune ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr); 1130d6ad1212SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1131d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1132d6ad1212SPeter Brune if (!fas->smoothd) { 11333de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1134d6ad1212SPeter Brune } 1135a3a80b83SMatthew G. Knepley *coarse = fas->smoothd; 1136d6ad1212SPeter Brune PetscFunctionReturn(0); 1137d6ad1212SPeter Brune } 1138928e959bSPeter Brune 1139928e959bSPeter Brune /*@ 1140928e959bSPeter Brune SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS 1141928e959bSPeter Brune 1142928e959bSPeter Brune Logically Collective on SNES 1143928e959bSPeter Brune 1144928e959bSPeter Brune Input Parameters: 1145928e959bSPeter Brune + snes - the multigrid context 1146928e959bSPeter Brune - swp - whether to downsweep or not 1147928e959bSPeter Brune 1148928e959bSPeter Brune Options Database Key: 1149928e959bSPeter Brune . -snes_fas_full_downsweep - Sets number of pre-smoothing steps 1150928e959bSPeter Brune 1151928e959bSPeter Brune Level: advanced 1152928e959bSPeter Brune 1153928e959bSPeter Brune .seealso: SNESFASSetNumberSmoothUp() 1154928e959bSPeter Brune @*/ 1155928e959bSPeter Brune PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp) 1156928e959bSPeter Brune { 1157f833ba53SLisandro Dalcin SNES_FAS *fas; 1158f833ba53SLisandro Dalcin PetscErrorCode ierr; 1159928e959bSPeter Brune 1160928e959bSPeter Brune PetscFunctionBegin; 1161f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 1162f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 1163928e959bSPeter Brune fas->full_downsweep = swp; 1164928e959bSPeter Brune if (fas->next) { 1165928e959bSPeter Brune ierr = SNESFASFullSetDownSweep(fas->next,swp);CHKERRQ(ierr); 1166928e959bSPeter Brune } 1167928e959bSPeter Brune PetscFunctionReturn(0); 1168928e959bSPeter Brune } 11696dfbafefSToby Isaac 11706dfbafefSToby Isaac /*@ 11716dfbafefSToby Isaac SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 11726dfbafefSToby Isaac 11736dfbafefSToby Isaac Logically Collective on SNES 11746dfbafefSToby Isaac 11756dfbafefSToby Isaac Input Parameters: 11766dfbafefSToby Isaac + snes - the multigrid context 11776dfbafefSToby Isaac - total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 11786dfbafefSToby Isaac 11796dfbafefSToby Isaac Options Database Key: 11806dfbafefSToby Isaac . -snes_fas_full_total 11816dfbafefSToby Isaac 11826dfbafefSToby Isaac Level: advanced 11836dfbafefSToby Isaac 11846dfbafefSToby Isaac Note: this option is only significant if the interpolation of a coarse correction (MatInterpolate()) is significantly different from total solution interpolation (DMInterpolateSolution()). 11856dfbafefSToby Isaac 11866dfbafefSToby Isaac .seealso: SNESFASSetNumberSmoothUp(), DMInterpolateSolution() 11876dfbafefSToby Isaac @*/ 11886dfbafefSToby Isaac PetscErrorCode SNESFASFullSetTotal(SNES snes,PetscBool total) 11896dfbafefSToby Isaac { 11906dfbafefSToby Isaac SNES_FAS *fas; 11916dfbafefSToby Isaac PetscErrorCode ierr; 11926dfbafefSToby Isaac 11936dfbafefSToby Isaac PetscFunctionBegin; 11946dfbafefSToby Isaac PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 11956dfbafefSToby Isaac fas = (SNES_FAS*)snes->data; 11966dfbafefSToby Isaac fas->full_total = total; 11976dfbafefSToby Isaac if (fas->next) { 11986dfbafefSToby Isaac ierr = SNESFASFullSetTotal(fas->next,total);CHKERRQ(ierr); 11996dfbafefSToby Isaac } 12006dfbafefSToby Isaac PetscFunctionReturn(0); 12016dfbafefSToby Isaac } 12026dfbafefSToby Isaac 12036dfbafefSToby Isaac /*@ 12046dfbafefSToby Isaac SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles 12056dfbafefSToby Isaac 12066dfbafefSToby Isaac Logically Collective on SNES 12076dfbafefSToby Isaac 12086dfbafefSToby Isaac Input Parameters: 12096dfbafefSToby Isaac . snes - the multigrid context 12106dfbafefSToby Isaac 12116dfbafefSToby Isaac Output: 12126dfbafefSToby Isaac . total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation) 12136dfbafefSToby Isaac 12146dfbafefSToby Isaac Options Database Key: 12156dfbafefSToby Isaac . -snes_fas_full_total 12166dfbafefSToby Isaac 12176dfbafefSToby Isaac Level: advanced 12186dfbafefSToby Isaac 12196dfbafefSToby Isaac .seealso: SNESFASSetNumberSmoothUp(), DMInterpolateSolution(), SNESFullSetTotal() 12206dfbafefSToby Isaac @*/ 12216dfbafefSToby Isaac PetscErrorCode SNESFASFullGetTotal(SNES snes,PetscBool *total) 12226dfbafefSToby Isaac { 12236dfbafefSToby Isaac SNES_FAS *fas; 12246dfbafefSToby Isaac 12256dfbafefSToby Isaac PetscFunctionBegin; 12266dfbafefSToby Isaac PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 12276dfbafefSToby Isaac fas = (SNES_FAS*)snes->data; 12286dfbafefSToby Isaac *total = fas->full_total; 12296dfbafefSToby Isaac PetscFunctionReturn(0); 12306dfbafefSToby Isaac } 1231