1a7b5fb5fSBarry Smith #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnes.h" I*/ 2ab8d36c9SPeter Brune 3ab8d36c9SPeter Brune 4ab8d36c9SPeter Brune extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES*); 5ab8d36c9SPeter Brune 6ab8d36c9SPeter Brune /* -------------- functions called on the fine level -------------- */ 7ab8d36c9SPeter Brune 8ab8d36c9SPeter Brune #undef __FUNCT__ 9ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetType" 10ab8d36c9SPeter Brune /*@ 11ab8d36c9SPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 12ab8d36c9SPeter Brune 13ab8d36c9SPeter Brune Logically Collective 14ab8d36c9SPeter Brune 15ab8d36c9SPeter Brune Input Parameters: 16583a1250SSatish Balay + snes - FAS context 1734d65b3cSPeter Brune - fastype - SNES_FAS_ADDITIVE, SNES_FAS_MULTIPLICATIVE, SNES_FAS_FULL, or SNES_FAS_KASKADE 18583a1250SSatish Balay 19583a1250SSatish Balay Level: intermediate 20ab8d36c9SPeter Brune 21ab8d36c9SPeter Brune .seealso: PCMGSetType() 22ab8d36c9SPeter Brune @*/ 23ab8d36c9SPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 24ab8d36c9SPeter Brune { 25ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 26ab8d36c9SPeter Brune PetscErrorCode ierr; 2722d28d08SBarry Smith 28ab8d36c9SPeter Brune PetscFunctionBegin; 29ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 30ab8d36c9SPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 31ab8d36c9SPeter Brune fas->fastype = fastype; 3222d28d08SBarry Smith if (fas->next) { 3322d28d08SBarry Smith ierr = SNESFASSetType(fas->next, fastype);CHKERRQ(ierr); 3422d28d08SBarry Smith } 35ab8d36c9SPeter Brune PetscFunctionReturn(0); 36ab8d36c9SPeter Brune } 37ab8d36c9SPeter Brune 38ab8d36c9SPeter Brune 39ab8d36c9SPeter Brune #undef __FUNCT__ 40ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetType" 41ab8d36c9SPeter Brune /*@ 42ab8d36c9SPeter Brune SNESFASGetType - Sets the update and correction type used for FAS. 43ab8d36c9SPeter Brune 44ab8d36c9SPeter Brune Logically Collective 45ab8d36c9SPeter Brune 46ab8d36c9SPeter Brune Input Parameters: 47ab8d36c9SPeter Brune . snes - FAS context 48ab8d36c9SPeter Brune 49ab8d36c9SPeter Brune Output Parameters: 50ab8d36c9SPeter Brune . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE 51ab8d36c9SPeter Brune 52583a1250SSatish Balay Level: intermediate 53583a1250SSatish Balay 54ab8d36c9SPeter Brune .seealso: PCMGSetType() 55ab8d36c9SPeter Brune @*/ 56ab8d36c9SPeter Brune PetscErrorCode SNESFASGetType(SNES snes,SNESFASType *fastype) 57ab8d36c9SPeter Brune { 58ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 59ab8d36c9SPeter Brune 60ab8d36c9SPeter Brune PetscFunctionBegin; 61ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 62ab8d36c9SPeter Brune PetscValidPointer(fastype, 2); 63ab8d36c9SPeter Brune *fastype = fas->fastype; 64ab8d36c9SPeter Brune PetscFunctionReturn(0); 65ab8d36c9SPeter Brune } 66ab8d36c9SPeter Brune 67ab8d36c9SPeter Brune #undef __FUNCT__ 68ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 69ab8d36c9SPeter Brune /*@C 70ab8d36c9SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 71ab8d36c9SPeter Brune Must be called before any other FAS routine. 72ab8d36c9SPeter Brune 73ab8d36c9SPeter Brune Input Parameters: 74ab8d36c9SPeter Brune + snes - the snes context 75ab8d36c9SPeter Brune . levels - the number of levels 76ab8d36c9SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 770298fd71SBarry Smith problems on smaller sets of processors. Use NULL_OBJECT for default in 78ab8d36c9SPeter Brune Fortran. 79ab8d36c9SPeter Brune 80ab8d36c9SPeter Brune Level: intermediate 81ab8d36c9SPeter Brune 82ab8d36c9SPeter Brune Notes: 83ab8d36c9SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 84ab8d36c9SPeter Brune for setting the level options rather than the -fas_coarse prefix. 85ab8d36c9SPeter Brune 86ab8d36c9SPeter Brune .keywords: FAS, MG, set, levels, multigrid 87ab8d36c9SPeter Brune 88ab8d36c9SPeter Brune .seealso: SNESFASGetLevels() 89ab8d36c9SPeter Brune @*/ 9022d28d08SBarry Smith PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) 9122d28d08SBarry Smith { 92ab8d36c9SPeter Brune PetscErrorCode ierr; 93ab8d36c9SPeter Brune PetscInt i; 94ab8d36c9SPeter Brune const char *optionsprefix; 95ab8d36c9SPeter Brune char tprefix[128]; 96ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 97ab8d36c9SPeter Brune SNES prevsnes; 98ab8d36c9SPeter Brune MPI_Comm comm; 9922d28d08SBarry Smith 100ab8d36c9SPeter Brune PetscFunctionBegin; 101ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 102ab8d36c9SPeter Brune if (levels == fas->levels) { 10322d28d08SBarry Smith if (!comms) PetscFunctionReturn(0); 104ab8d36c9SPeter Brune } 105ab8d36c9SPeter Brune /* user has changed the number of levels; reset */ 106ab8d36c9SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 107ab8d36c9SPeter Brune /* destroy any coarser levels if necessary */ 108ab8d36c9SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 1090298fd71SBarry Smith fas->next = NULL; 1100298fd71SBarry Smith fas->previous = NULL; 111ab8d36c9SPeter Brune prevsnes = snes; 112ab8d36c9SPeter Brune /* setup the finest level */ 113ab8d36c9SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 114ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1);CHKERRQ(ierr); 115ab8d36c9SPeter Brune for (i = levels - 1; i >= 0; i--) { 116ab8d36c9SPeter Brune if (comms) comm = comms[i]; 117ab8d36c9SPeter Brune fas->level = i; 118ab8d36c9SPeter Brune fas->levels = levels; 119ab8d36c9SPeter Brune fas->fine = snes; 1200298fd71SBarry Smith fas->next = NULL; 121ab8d36c9SPeter Brune if (i > 0) { 122ab8d36c9SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 123e964f0dbSPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 124ab8d36c9SPeter Brune sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level); 125ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(fas->next,optionsprefix);CHKERRQ(ierr); 126e964f0dbSPeter Brune ierr = SNESAppendOptionsPrefix(fas->next,tprefix);CHKERRQ(ierr); 127ab8d36c9SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 128ab8d36c9SPeter Brune ierr = SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);CHKERRQ(ierr); 129ab8d36c9SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 130ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1);CHKERRQ(ierr); 1311aa26658SKarl Rupp 132ab8d36c9SPeter Brune ((SNES_FAS*)fas->next->data)->previous = prevsnes; 1331aa26658SKarl Rupp 134ab8d36c9SPeter Brune prevsnes = fas->next; 135ab8d36c9SPeter Brune fas = (SNES_FAS*)prevsnes->data; 136ab8d36c9SPeter Brune } 137ab8d36c9SPeter Brune } 138ab8d36c9SPeter Brune PetscFunctionReturn(0); 139ab8d36c9SPeter Brune } 140ab8d36c9SPeter Brune 141ab8d36c9SPeter Brune 142ab8d36c9SPeter Brune #undef __FUNCT__ 143ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 144ab8d36c9SPeter Brune /*@ 145ab8d36c9SPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 146ab8d36c9SPeter Brune 147ab8d36c9SPeter Brune Input Parameter: 148ab8d36c9SPeter Brune . snes - the nonlinear solver context 149ab8d36c9SPeter Brune 150ab8d36c9SPeter Brune Output parameter: 151ab8d36c9SPeter Brune . levels - the number of levels 152ab8d36c9SPeter Brune 153ab8d36c9SPeter Brune Level: advanced 154ab8d36c9SPeter Brune 155ab8d36c9SPeter Brune .keywords: MG, get, levels, multigrid 156ab8d36c9SPeter Brune 157ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 158ab8d36c9SPeter Brune @*/ 15922d28d08SBarry Smith PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels) 16022d28d08SBarry Smith { 161ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 1625fd66863SKarl Rupp 163ab8d36c9SPeter Brune PetscFunctionBegin; 164ab8d36c9SPeter Brune *levels = fas->levels; 165ab8d36c9SPeter Brune PetscFunctionReturn(0); 166ab8d36c9SPeter Brune } 167ab8d36c9SPeter Brune 168ab8d36c9SPeter Brune 169ab8d36c9SPeter Brune #undef __FUNCT__ 170ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetCycleSNES" 171ab8d36c9SPeter Brune /*@ 172ab8d36c9SPeter Brune SNESFASGetCycleSNES - Gets the SNES corresponding to a particular 173ab8d36c9SPeter Brune level of the FAS hierarchy. 174ab8d36c9SPeter Brune 175ab8d36c9SPeter Brune Input Parameters: 176ab8d36c9SPeter Brune + snes - the multigrid context 177ab8d36c9SPeter Brune level - the level to get 178ab8d36c9SPeter Brune - lsnes - whether to use the nonlinear smoother or not 179ab8d36c9SPeter Brune 180ab8d36c9SPeter Brune Level: advanced 181ab8d36c9SPeter Brune 182ab8d36c9SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 183ab8d36c9SPeter Brune 184ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 185ab8d36c9SPeter Brune @*/ 18622d28d08SBarry Smith PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes) 18722d28d08SBarry Smith { 188ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 189ab8d36c9SPeter Brune PetscInt i; 190ab8d36c9SPeter Brune 191ab8d36c9SPeter Brune PetscFunctionBegin; 192ce94432eSBarry Smith if (level > fas->levels-1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels); 193ce94432eSBarry Smith if (fas->level != fas->levels - 1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level); 194ab8d36c9SPeter Brune 195ab8d36c9SPeter Brune *lsnes = snes; 196ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) { 197ab8d36c9SPeter Brune *lsnes = fas->next; 198ab8d36c9SPeter Brune fas = (SNES_FAS*)(*lsnes)->data; 199ab8d36c9SPeter Brune } 200ce94432eSBarry Smith if (fas->level != level) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 201ab8d36c9SPeter Brune PetscFunctionReturn(0); 202ab8d36c9SPeter Brune } 203ab8d36c9SPeter Brune 204ab8d36c9SPeter Brune #undef __FUNCT__ 205ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 206ab8d36c9SPeter Brune /*@ 207ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 208ab8d36c9SPeter Brune use on all levels. 209ab8d36c9SPeter Brune 210ab8d36c9SPeter Brune Logically Collective on SNES 211ab8d36c9SPeter Brune 212ab8d36c9SPeter Brune Input Parameters: 213ab8d36c9SPeter Brune + snes - the multigrid context 214ab8d36c9SPeter Brune - n - the number of smoothing steps 215ab8d36c9SPeter Brune 216ab8d36c9SPeter Brune Options Database Key: 217ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 218ab8d36c9SPeter Brune 219ab8d36c9SPeter Brune Level: advanced 220ab8d36c9SPeter Brune 221ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 222ab8d36c9SPeter Brune 223ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 224ab8d36c9SPeter Brune @*/ 22522d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) 22622d28d08SBarry Smith { 227ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 22822d28d08SBarry Smith PetscErrorCode ierr; 22922d28d08SBarry Smith 230ab8d36c9SPeter Brune PetscFunctionBegin; 231ab8d36c9SPeter Brune fas->max_up_it = n; 232656ede7eSPeter Brune if (!fas->smoothu && fas->level != 0) { 23322d28d08SBarry Smith ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 234ab8d36c9SPeter Brune } 23522d28d08SBarry Smith if (fas->smoothu) { 23622d28d08SBarry Smith ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr); 23722d28d08SBarry Smith } 238ab8d36c9SPeter Brune if (fas->next) { 239ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 240ab8d36c9SPeter Brune } 241ab8d36c9SPeter Brune PetscFunctionReturn(0); 242ab8d36c9SPeter Brune } 243ab8d36c9SPeter Brune 244ab8d36c9SPeter Brune #undef __FUNCT__ 245ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 246ab8d36c9SPeter Brune /*@ 247ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 248ab8d36c9SPeter Brune use on all levels. 249ab8d36c9SPeter Brune 250ab8d36c9SPeter Brune Logically Collective on SNES 251ab8d36c9SPeter Brune 252ab8d36c9SPeter Brune Input Parameters: 253ab8d36c9SPeter Brune + snes - the multigrid context 254ab8d36c9SPeter Brune - n - the number of smoothing steps 255ab8d36c9SPeter Brune 256ab8d36c9SPeter Brune Options Database Key: 257ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 258ab8d36c9SPeter Brune 259ab8d36c9SPeter Brune Level: advanced 260ab8d36c9SPeter Brune 261ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 262ab8d36c9SPeter Brune 263ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 264ab8d36c9SPeter Brune @*/ 26522d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) 26622d28d08SBarry Smith { 267ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 268ab8d36c9SPeter Brune PetscErrorCode ierr = 0; 26922d28d08SBarry Smith 270ab8d36c9SPeter Brune PetscFunctionBegin; 271ab8d36c9SPeter Brune if (!fas->smoothd) { 27222d28d08SBarry Smith ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 273ab8d36c9SPeter Brune } 274ab8d36c9SPeter Brune ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr); 2751aa26658SKarl Rupp 276ab8d36c9SPeter Brune fas->max_down_it = n; 277ab8d36c9SPeter Brune if (fas->next) { 278ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 279ab8d36c9SPeter Brune } 280ab8d36c9SPeter Brune PetscFunctionReturn(0); 281ab8d36c9SPeter Brune } 282ab8d36c9SPeter Brune 283ab8d36c9SPeter Brune 284ab8d36c9SPeter Brune #undef __FUNCT__ 28587f44e3fSPeter Brune #define __FUNCT__ "SNESFASSetContinuation" 28687f44e3fSPeter Brune /*@ 28787f44e3fSPeter Brune SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep 28887f44e3fSPeter Brune 28987f44e3fSPeter Brune Logically Collective on SNES 29087f44e3fSPeter Brune 29187f44e3fSPeter Brune Input Parameters: 29287f44e3fSPeter Brune + snes - the multigrid context 29387f44e3fSPeter Brune - n - the number of smoothing steps 29487f44e3fSPeter Brune 29587f44e3fSPeter Brune Options Database Key: 29687f44e3fSPeter Brune . -snes_fas_continuation - sets continuation to true 29787f44e3fSPeter Brune 29887f44e3fSPeter Brune Level: advanced 29987f44e3fSPeter Brune 30087f44e3fSPeter Brune Notes: This sets the prefix on the upsweep smoothers to -fas_continuation 30187f44e3fSPeter Brune 30287f44e3fSPeter Brune .keywords: FAS, MG, smoother, continuation 30387f44e3fSPeter Brune 30487f44e3fSPeter Brune .seealso: SNESFAS 30587f44e3fSPeter Brune @*/ 30687f44e3fSPeter Brune PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation) 30787f44e3fSPeter Brune { 30887f44e3fSPeter Brune const char *optionsprefix; 30987f44e3fSPeter Brune char tprefix[128]; 31087f44e3fSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 31187f44e3fSPeter Brune PetscErrorCode ierr = 0; 31287f44e3fSPeter Brune 31387f44e3fSPeter Brune PetscFunctionBegin; 31487f44e3fSPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 31587f44e3fSPeter Brune if (!fas->smoothu) { 31687f44e3fSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 31787f44e3fSPeter Brune } 31887f44e3fSPeter Brune sprintf(tprefix,"fas_levels_continuation_"); 31987f44e3fSPeter Brune ierr = SNESSetOptionsPrefix(fas->smoothu, optionsprefix);CHKERRQ(ierr); 32087f44e3fSPeter Brune ierr = SNESAppendOptionsPrefix(fas->smoothu, tprefix);CHKERRQ(ierr); 32187f44e3fSPeter Brune ierr = SNESSetType(fas->smoothu,SNESNEWTONLS);CHKERRQ(ierr); 32287f44e3fSPeter Brune ierr = SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100);CHKERRQ(ierr); 32387f44e3fSPeter Brune fas->continuation = continuation; 32487f44e3fSPeter Brune if (fas->next) { 32587f44e3fSPeter Brune ierr = SNESFASSetContinuation(fas->next,continuation);CHKERRQ(ierr); 32687f44e3fSPeter Brune } 32787f44e3fSPeter Brune PetscFunctionReturn(0); 32887f44e3fSPeter Brune } 32987f44e3fSPeter Brune 33087f44e3fSPeter Brune 33187f44e3fSPeter Brune #undef __FUNCT__ 332ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetCycles" 333ab8d36c9SPeter Brune /*@ 334ab8d36c9SPeter Brune SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 335ab8d36c9SPeter Brune complicated cycling. 336ab8d36c9SPeter Brune 337ab8d36c9SPeter Brune Logically Collective on SNES 338ab8d36c9SPeter Brune 339ab8d36c9SPeter Brune Input Parameters: 340ab8d36c9SPeter Brune + snes - the multigrid context 341ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 342ab8d36c9SPeter Brune 343ab8d36c9SPeter Brune Options Database Key: 344e1bc860dSBarry Smith . -snes_fas_cycles 1 or 2 345ab8d36c9SPeter Brune 346ab8d36c9SPeter Brune Level: advanced 347ab8d36c9SPeter Brune 348ab8d36c9SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 349ab8d36c9SPeter Brune 350ab8d36c9SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 351ab8d36c9SPeter Brune @*/ 35222d28d08SBarry Smith PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) 35322d28d08SBarry Smith { 354ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 355ab8d36c9SPeter Brune PetscErrorCode ierr; 356ab8d36c9SPeter Brune PetscBool isFine; 35722d28d08SBarry Smith 358ab8d36c9SPeter Brune PetscFunctionBegin; 35922d28d08SBarry Smith ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 3601aa26658SKarl Rupp 361ab8d36c9SPeter Brune fas->n_cycles = cycles; 3621aa26658SKarl Rupp if (!isFine) { 363ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 3641aa26658SKarl Rupp } 365ab8d36c9SPeter Brune if (fas->next) { 366ab8d36c9SPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 367ab8d36c9SPeter Brune } 368ab8d36c9SPeter Brune PetscFunctionReturn(0); 369ab8d36c9SPeter Brune } 370ab8d36c9SPeter Brune 371c8c899caSPeter Brune 372c8c899caSPeter Brune #undef __FUNCT__ 373c8c899caSPeter Brune #define __FUNCT__ "SNESFASSetMonitor" 374c8c899caSPeter Brune /*@ 375c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring 376c8c899caSPeter Brune 377c8c899caSPeter Brune Logically Collective on SNES 378c8c899caSPeter Brune 379c8c899caSPeter Brune Input Parameters: 380c8c899caSPeter Brune + snes - the FAS context 381*d142ab34SLawrence Mitchell . vf - viewer and format structure (may be NULL if flg is FALSE) 382c8c899caSPeter Brune - flg - monitor or not 383c8c899caSPeter Brune 384c8c899caSPeter Brune Level: advanced 385c8c899caSPeter Brune 386c8c899caSPeter Brune .keywords: FAS, monitor 387c8c899caSPeter Brune 388c8c899caSPeter Brune .seealso: SNESFASSetCyclesOnLevel() 389c8c899caSPeter Brune @*/ 390*d142ab34SLawrence Mitchell PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg) 39122d28d08SBarry Smith { 392c8c899caSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 393c8c899caSPeter Brune PetscErrorCode ierr; 394c8c899caSPeter Brune PetscBool isFine; 395c8c899caSPeter Brune PetscInt i, levels = fas->levels; 396c8c899caSPeter Brune SNES levelsnes; 39722d28d08SBarry Smith 398c8c899caSPeter Brune PetscFunctionBegin; 399c8c899caSPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 400c8c899caSPeter Brune if (isFine) { 401c8c899caSPeter Brune for (i = 0; i < levels; i++) { 40222d28d08SBarry Smith ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 403c8c899caSPeter Brune fas = (SNES_FAS*)levelsnes->data; 404c8c899caSPeter Brune if (flg) { 405c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */ 406c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 407*d142ab34SLawrence Mitchell /* Only register destroy on finest level */ 408*d142ab34SLawrence Mitchell ierr = SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL));CHKERRQ(ierr); 4091aa26658SKarl Rupp } else if (i != fas->levels - 1) { 410c8c899caSPeter Brune /* unset the monitors on the coarse levels */ 411c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 412c8c899caSPeter Brune } 413c8c899caSPeter Brune } 414c8c899caSPeter Brune } 415c8c899caSPeter Brune PetscFunctionReturn(0); 416c8c899caSPeter Brune } 417c8c899caSPeter Brune 418ab8d36c9SPeter Brune #undef __FUNCT__ 4190dd27c6cSPeter Brune #define __FUNCT__ "SNESFASSetLog" 4200dd27c6cSPeter Brune /*@ 4210dd27c6cSPeter Brune SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels 4220dd27c6cSPeter Brune 4230dd27c6cSPeter Brune Logically Collective on SNES 4240dd27c6cSPeter Brune 4250dd27c6cSPeter Brune Input Parameters: 4260dd27c6cSPeter Brune + snes - the FAS context 4270dd27c6cSPeter Brune - flg - monitor or not 4280dd27c6cSPeter Brune 4290dd27c6cSPeter Brune Level: advanced 4300dd27c6cSPeter Brune 4310dd27c6cSPeter Brune .keywords: FAS, logging 4320dd27c6cSPeter Brune 4330dd27c6cSPeter Brune .seealso: SNESFASSetMonitor() 4340dd27c6cSPeter Brune @*/ 4350dd27c6cSPeter Brune PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) 4360dd27c6cSPeter Brune { 4370dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 4380dd27c6cSPeter Brune PetscErrorCode ierr; 4390dd27c6cSPeter Brune PetscBool isFine; 4400dd27c6cSPeter Brune PetscInt i, levels = fas->levels; 4410dd27c6cSPeter Brune SNES levelsnes; 4420dd27c6cSPeter Brune char eventname[128]; 4430dd27c6cSPeter Brune 4440dd27c6cSPeter Brune PetscFunctionBegin; 4450dd27c6cSPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 4460dd27c6cSPeter Brune if (isFine) { 4470dd27c6cSPeter Brune for (i = 0; i < levels; i++) { 4480dd27c6cSPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 4490dd27c6cSPeter Brune fas = (SNES_FAS*)levelsnes->data; 4500dd27c6cSPeter Brune if (flg) { 4510dd27c6cSPeter Brune sprintf(eventname,"FASSetup %d",(int)i); 4520dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr); 4530dd27c6cSPeter Brune sprintf(eventname,"FASSmooth %d",(int)i); 4540dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr); 4550dd27c6cSPeter Brune sprintf(eventname,"FASResid %d",(int)i); 4560dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr); 4570dd27c6cSPeter Brune if (i) { 4580dd27c6cSPeter Brune sprintf(eventname,"FASInterp %d",(int)i); 4590dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr); 4600dd27c6cSPeter Brune } 4610dd27c6cSPeter Brune } else { 4620298fd71SBarry Smith fas->eventsmoothsetup = 0; 4630298fd71SBarry Smith fas->eventsmoothsolve = 0; 4640298fd71SBarry Smith fas->eventresidual = 0; 4650298fd71SBarry Smith fas->eventinterprestrict = 0; 4660dd27c6cSPeter Brune } 4670dd27c6cSPeter Brune } 4680dd27c6cSPeter Brune } 4690dd27c6cSPeter Brune PetscFunctionReturn(0); 4700dd27c6cSPeter Brune } 4710dd27c6cSPeter Brune 4720dd27c6cSPeter Brune #undef __FUNCT__ 473ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleCreateSmoother_Private" 474ab8d36c9SPeter Brune /* 475ab8d36c9SPeter Brune Creates the default smoother type. 476ab8d36c9SPeter Brune 47704d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level. 478ab8d36c9SPeter Brune 479ab8d36c9SPeter Brune */ 48022d28d08SBarry Smith PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) 48122d28d08SBarry Smith { 482ab8d36c9SPeter Brune SNES_FAS *fas; 483ab8d36c9SPeter Brune const char *optionsprefix; 484ab8d36c9SPeter Brune char tprefix[128]; 485ab8d36c9SPeter Brune PetscErrorCode ierr; 486ab8d36c9SPeter Brune SNES nsmooth; 48722d28d08SBarry Smith 488ab8d36c9SPeter Brune PetscFunctionBegin; 489ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 490ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 491ab8d36c9SPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 492ab8d36c9SPeter Brune /* create the default smoother */ 493ce94432eSBarry Smith ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);CHKERRQ(ierr); 494ab8d36c9SPeter Brune if (fas->level == 0) { 495ab8d36c9SPeter Brune sprintf(tprefix,"fas_coarse_"); 496ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 497ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 49804d7464bSBarry Smith ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr); 499e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr); 500ab8d36c9SPeter Brune } else { 501ab8d36c9SPeter Brune sprintf(tprefix,"fas_levels_%d_",(int)fas->level); 502ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 503ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 504ab8d36c9SPeter Brune ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr); 505e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr); 506ab8d36c9SPeter Brune } 507ab8d36c9SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 5083bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);CHKERRQ(ierr); 509f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr); 510ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level);CHKERRQ(ierr); 511ab8d36c9SPeter Brune *smooth = nsmooth; 512ab8d36c9SPeter Brune PetscFunctionReturn(0); 513ab8d36c9SPeter Brune } 514ab8d36c9SPeter Brune 515ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */ 516ab8d36c9SPeter Brune 517ab8d36c9SPeter Brune #undef __FUNCT__ 518ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleSetCycles" 519ab8d36c9SPeter Brune /*@ 520ab8d36c9SPeter Brune SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 521ab8d36c9SPeter Brune 522ab8d36c9SPeter Brune Logically Collective on SNES 523ab8d36c9SPeter Brune 524ab8d36c9SPeter Brune Input Parameters: 525ab8d36c9SPeter Brune + snes - the multigrid context 526ab8d36c9SPeter Brune . level - the level to set the number of cycles on 527ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 528ab8d36c9SPeter Brune 529ab8d36c9SPeter Brune Level: advanced 530ab8d36c9SPeter Brune 531ab8d36c9SPeter Brune .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid 532ab8d36c9SPeter Brune 533ab8d36c9SPeter Brune .seealso: SNESFASSetCycles() 534ab8d36c9SPeter Brune @*/ 53522d28d08SBarry Smith PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) 53622d28d08SBarry Smith { 537ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 538ab8d36c9SPeter Brune PetscErrorCode ierr; 53922d28d08SBarry Smith 540ab8d36c9SPeter Brune PetscFunctionBegin; 541ab8d36c9SPeter Brune fas->n_cycles = cycles; 542ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 543ab8d36c9SPeter Brune PetscFunctionReturn(0); 544ab8d36c9SPeter Brune } 545ab8d36c9SPeter Brune 546ab8d36c9SPeter Brune 547ab8d36c9SPeter Brune #undef __FUNCT__ 548ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmoother" 549ab8d36c9SPeter Brune /*@ 550ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 551ab8d36c9SPeter Brune 552ab8d36c9SPeter Brune Logically Collective on SNES 553ab8d36c9SPeter Brune 554ab8d36c9SPeter Brune Input Parameters: 555ab8d36c9SPeter Brune . snes - the multigrid context 556ab8d36c9SPeter Brune 557ab8d36c9SPeter Brune Output Parameters: 558ab8d36c9SPeter Brune . smooth - the smoother 559ab8d36c9SPeter Brune 560ab8d36c9SPeter Brune Level: advanced 561ab8d36c9SPeter Brune 562ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 563ab8d36c9SPeter Brune 564ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown() 565ab8d36c9SPeter Brune @*/ 566ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 567ab8d36c9SPeter Brune { 568ab8d36c9SPeter Brune SNES_FAS *fas; 5695fd66863SKarl Rupp 570ab8d36c9SPeter Brune PetscFunctionBegin; 571ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 572ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 573ab8d36c9SPeter Brune *smooth = fas->smoothd; 574ab8d36c9SPeter Brune PetscFunctionReturn(0); 575ab8d36c9SPeter Brune } 576ab8d36c9SPeter Brune #undef __FUNCT__ 577ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherUp" 578ab8d36c9SPeter Brune /*@ 579ab8d36c9SPeter Brune SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 580ab8d36c9SPeter Brune 581ab8d36c9SPeter Brune Logically Collective on SNES 582ab8d36c9SPeter Brune 583ab8d36c9SPeter Brune Input Parameters: 584ab8d36c9SPeter Brune . snes - the multigrid context 585ab8d36c9SPeter Brune 586ab8d36c9SPeter Brune Output Parameters: 587ab8d36c9SPeter Brune . smoothu - the smoother 588ab8d36c9SPeter Brune 589ab8d36c9SPeter Brune Notes: 590ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent 591ab8d36c9SPeter Brune default behavior in the process of the solve. 592ab8d36c9SPeter Brune 593ab8d36c9SPeter Brune Level: advanced 594ab8d36c9SPeter Brune 595ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 596ab8d36c9SPeter Brune 597ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown() 598ab8d36c9SPeter Brune @*/ 599ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 600ab8d36c9SPeter Brune { 601ab8d36c9SPeter Brune SNES_FAS *fas; 6025fd66863SKarl Rupp 603ab8d36c9SPeter Brune PetscFunctionBegin; 604ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 605ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 6061aa26658SKarl Rupp if (!fas->smoothu) *smoothu = fas->smoothd; 6071aa26658SKarl Rupp else *smoothu = fas->smoothu; 608ab8d36c9SPeter Brune PetscFunctionReturn(0); 609ab8d36c9SPeter Brune } 610ab8d36c9SPeter Brune 611ab8d36c9SPeter Brune #undef __FUNCT__ 612ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherDown" 613ab8d36c9SPeter Brune /*@ 614ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 615ab8d36c9SPeter Brune 616ab8d36c9SPeter Brune Logically Collective on SNES 617ab8d36c9SPeter Brune 618ab8d36c9SPeter Brune Input Parameters: 619ab8d36c9SPeter Brune . snes - the multigrid context 620ab8d36c9SPeter Brune 621ab8d36c9SPeter Brune Output Parameters: 622ab8d36c9SPeter Brune . smoothd - the smoother 623ab8d36c9SPeter Brune 624ab8d36c9SPeter Brune Level: advanced 625ab8d36c9SPeter Brune 626ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 627ab8d36c9SPeter Brune 628ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 629ab8d36c9SPeter Brune @*/ 630ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 631ab8d36c9SPeter Brune { 632ab8d36c9SPeter Brune SNES_FAS *fas; 6335fd66863SKarl Rupp 634ab8d36c9SPeter Brune PetscFunctionBegin; 635ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 636ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 637ab8d36c9SPeter Brune *smoothd = fas->smoothd; 638ab8d36c9SPeter Brune PetscFunctionReturn(0); 639ab8d36c9SPeter Brune } 640ab8d36c9SPeter Brune 641ab8d36c9SPeter Brune 642ab8d36c9SPeter Brune #undef __FUNCT__ 643ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetCorrection" 644ab8d36c9SPeter Brune /*@ 645ab8d36c9SPeter Brune SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 646ab8d36c9SPeter Brune 647ab8d36c9SPeter Brune Logically Collective on SNES 648ab8d36c9SPeter Brune 649ab8d36c9SPeter Brune Input Parameters: 650ab8d36c9SPeter Brune . snes - the multigrid context 651ab8d36c9SPeter Brune 652ab8d36c9SPeter Brune Output Parameters: 653ab8d36c9SPeter Brune . correction - the coarse correction on this level 654ab8d36c9SPeter Brune 655ab8d36c9SPeter Brune Notes: 6560298fd71SBarry Smith Returns NULL on the coarsest level. 657ab8d36c9SPeter Brune 658ab8d36c9SPeter Brune Level: advanced 659ab8d36c9SPeter Brune 660ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 661ab8d36c9SPeter Brune 662ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 663ab8d36c9SPeter Brune @*/ 664ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 665ab8d36c9SPeter Brune { 666ab8d36c9SPeter Brune SNES_FAS *fas; 6675fd66863SKarl Rupp 668ab8d36c9SPeter Brune PetscFunctionBegin; 669ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 670ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 671ab8d36c9SPeter Brune *correction = fas->next; 672ab8d36c9SPeter Brune PetscFunctionReturn(0); 673ab8d36c9SPeter Brune } 674ab8d36c9SPeter Brune 675ab8d36c9SPeter Brune #undef __FUNCT__ 676ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInterpolation" 677ab8d36c9SPeter Brune /*@ 678ab8d36c9SPeter Brune SNESFASCycleGetInterpolation - Gets the interpolation on this level 679ab8d36c9SPeter Brune 680ab8d36c9SPeter Brune Logically Collective on SNES 681ab8d36c9SPeter Brune 682ab8d36c9SPeter Brune Input Parameters: 683ab8d36c9SPeter Brune . snes - the multigrid context 684ab8d36c9SPeter Brune 685ab8d36c9SPeter Brune Output Parameters: 686ab8d36c9SPeter Brune . mat - the interpolation operator on this level 687ab8d36c9SPeter Brune 688ab8d36c9SPeter Brune Level: developer 689ab8d36c9SPeter Brune 690ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 691ab8d36c9SPeter Brune 692ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 693ab8d36c9SPeter Brune @*/ 694ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 695ab8d36c9SPeter Brune { 696ab8d36c9SPeter Brune SNES_FAS *fas; 6975fd66863SKarl Rupp 698ab8d36c9SPeter Brune PetscFunctionBegin; 699ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 700ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 701ab8d36c9SPeter Brune *mat = fas->interpolate; 702ab8d36c9SPeter Brune PetscFunctionReturn(0); 703ab8d36c9SPeter Brune } 704ab8d36c9SPeter Brune 705ab8d36c9SPeter Brune 706ab8d36c9SPeter Brune #undef __FUNCT__ 707ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRestriction" 708ab8d36c9SPeter Brune /*@ 709ab8d36c9SPeter Brune SNESFASCycleGetRestriction - Gets the restriction on this level 710ab8d36c9SPeter Brune 711ab8d36c9SPeter Brune Logically Collective on SNES 712ab8d36c9SPeter Brune 713ab8d36c9SPeter Brune Input Parameters: 714ab8d36c9SPeter Brune . snes - the multigrid context 715ab8d36c9SPeter Brune 716ab8d36c9SPeter Brune Output Parameters: 717ab8d36c9SPeter Brune . mat - the restriction operator on this level 718ab8d36c9SPeter Brune 719ab8d36c9SPeter Brune Level: developer 720ab8d36c9SPeter Brune 721ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 722ab8d36c9SPeter Brune 723ab8d36c9SPeter Brune .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation() 724ab8d36c9SPeter Brune @*/ 725ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 726ab8d36c9SPeter Brune { 727ab8d36c9SPeter Brune SNES_FAS *fas; 7285fd66863SKarl Rupp 729ab8d36c9SPeter Brune PetscFunctionBegin; 730ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 731ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 732ab8d36c9SPeter Brune *mat = fas->restrct; 733ab8d36c9SPeter Brune PetscFunctionReturn(0); 734ab8d36c9SPeter Brune } 735ab8d36c9SPeter Brune 736ab8d36c9SPeter Brune 737ab8d36c9SPeter Brune #undef __FUNCT__ 738ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInjection" 739ab8d36c9SPeter Brune /*@ 740ab8d36c9SPeter Brune SNESFASCycleGetInjection - Gets the injection on this level 741ab8d36c9SPeter Brune 742ab8d36c9SPeter Brune Logically Collective on SNES 743ab8d36c9SPeter Brune 744ab8d36c9SPeter Brune Input Parameters: 745ab8d36c9SPeter Brune . snes - the multigrid context 746ab8d36c9SPeter Brune 747ab8d36c9SPeter Brune Output Parameters: 748ab8d36c9SPeter Brune . mat - the restriction operator on this level 749ab8d36c9SPeter Brune 750ab8d36c9SPeter Brune Level: developer 751ab8d36c9SPeter Brune 752ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 753ab8d36c9SPeter Brune 754ab8d36c9SPeter Brune .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction() 755ab8d36c9SPeter Brune @*/ 756ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 757ab8d36c9SPeter Brune { 758ab8d36c9SPeter Brune SNES_FAS *fas; 7595fd66863SKarl Rupp 760ab8d36c9SPeter Brune PetscFunctionBegin; 761ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 762ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 763ab8d36c9SPeter Brune *mat = fas->inject; 764ab8d36c9SPeter Brune PetscFunctionReturn(0); 765ab8d36c9SPeter Brune } 766ab8d36c9SPeter Brune 767ab8d36c9SPeter Brune #undef __FUNCT__ 768ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRScale" 769ab8d36c9SPeter Brune /*@ 770ab8d36c9SPeter Brune SNESFASCycleGetRScale - Gets the injection on this level 771ab8d36c9SPeter Brune 772ab8d36c9SPeter Brune Logically Collective on SNES 773ab8d36c9SPeter Brune 774ab8d36c9SPeter Brune Input Parameters: 775ab8d36c9SPeter Brune . snes - the multigrid context 776ab8d36c9SPeter Brune 777ab8d36c9SPeter Brune Output Parameters: 778ab8d36c9SPeter Brune . mat - the restriction operator on this level 779ab8d36c9SPeter Brune 780ab8d36c9SPeter Brune Level: developer 781ab8d36c9SPeter Brune 782ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 783ab8d36c9SPeter Brune 784ab8d36c9SPeter Brune .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale() 785ab8d36c9SPeter Brune @*/ 786ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 787ab8d36c9SPeter Brune { 788ab8d36c9SPeter Brune SNES_FAS *fas; 7895fd66863SKarl Rupp 790ab8d36c9SPeter Brune PetscFunctionBegin; 791ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 792ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 793ab8d36c9SPeter Brune *vec = fas->rscale; 794ab8d36c9SPeter Brune PetscFunctionReturn(0); 795ab8d36c9SPeter Brune } 796ab8d36c9SPeter Brune 797ab8d36c9SPeter Brune #undef __FUNCT__ 798ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleIsFine" 799ab8d36c9SPeter Brune /*@ 800ab8d36c9SPeter Brune SNESFASCycleIsFine - Determines if a given cycle is the fine level. 801ab8d36c9SPeter Brune 802ab8d36c9SPeter Brune Logically Collective on SNES 803ab8d36c9SPeter Brune 804ab8d36c9SPeter Brune Input Parameters: 805ab8d36c9SPeter Brune . snes - the FAS context 806ab8d36c9SPeter Brune 807ab8d36c9SPeter Brune Output Parameters: 808ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not 809ab8d36c9SPeter Brune 810ab8d36c9SPeter Brune Level: advanced 811ab8d36c9SPeter Brune 812ab8d36c9SPeter Brune .keywords: SNES, FAS 813ab8d36c9SPeter Brune 814ab8d36c9SPeter Brune .seealso: SNESFASSetLevels() 815ab8d36c9SPeter Brune @*/ 816ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 817ab8d36c9SPeter Brune { 818ab8d36c9SPeter Brune SNES_FAS *fas; 8195fd66863SKarl Rupp 820ab8d36c9SPeter Brune PetscFunctionBegin; 821ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 822ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 8231aa26658SKarl Rupp if (fas->level == fas->levels - 1) *flg = PETSC_TRUE; 8241aa26658SKarl Rupp else *flg = PETSC_FALSE; 825ab8d36c9SPeter Brune PetscFunctionReturn(0); 826ab8d36c9SPeter Brune } 827ab8d36c9SPeter Brune 828ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */ 829ab8d36c9SPeter Brune 830ab8d36c9SPeter Brune #undef __FUNCT__ 831ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 832ab8d36c9SPeter Brune /*@ 833ab8d36c9SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 834ab8d36c9SPeter Brune interpolation from l-1 to the lth level 835ab8d36c9SPeter Brune 836ab8d36c9SPeter Brune Input Parameters: 837ab8d36c9SPeter Brune + snes - the multigrid context 838ab8d36c9SPeter Brune . mat - the interpolation operator 839ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 840ab8d36c9SPeter Brune 841ab8d36c9SPeter Brune Level: advanced 842ab8d36c9SPeter Brune 843ab8d36c9SPeter Brune Notes: 844ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction 845ab8d36c9SPeter Brune for the same level. 846ab8d36c9SPeter Brune 847ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 848ab8d36c9SPeter Brune out from the matrix size which one it is. 849ab8d36c9SPeter Brune 850ab8d36c9SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 851ab8d36c9SPeter Brune 852ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 853ab8d36c9SPeter Brune @*/ 8540adebc6cSBarry Smith PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) 8550adebc6cSBarry Smith { 85622d28d08SBarry Smith SNES_FAS *fas; 857ab8d36c9SPeter Brune PetscErrorCode ierr; 858ab8d36c9SPeter Brune SNES levelsnes; 85922d28d08SBarry Smith 860ab8d36c9SPeter Brune PetscFunctionBegin; 861ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 862ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 863ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 864ab8d36c9SPeter Brune ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 8651aa26658SKarl Rupp 866ab8d36c9SPeter Brune fas->interpolate = mat; 867ab8d36c9SPeter Brune PetscFunctionReturn(0); 868ab8d36c9SPeter Brune } 869ab8d36c9SPeter Brune 870ab8d36c9SPeter Brune #undef __FUNCT__ 871ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInterpolation" 872ab8d36c9SPeter Brune /*@ 873ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the 874ab8d36c9SPeter Brune interpolation from l-1 to the lth level 875ab8d36c9SPeter Brune 876ab8d36c9SPeter Brune Input Parameters: 877ab8d36c9SPeter Brune + snes - the multigrid context 878ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 879ab8d36c9SPeter Brune 880ab8d36c9SPeter Brune Output Parameters: 881ab8d36c9SPeter Brune . mat - the interpolation operator 882ab8d36c9SPeter Brune 883ab8d36c9SPeter Brune Level: advanced 884ab8d36c9SPeter Brune 885ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, interpolate, level 886ab8d36c9SPeter Brune 887ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale() 888ab8d36c9SPeter Brune @*/ 88922d28d08SBarry Smith PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) 89022d28d08SBarry Smith { 89122d28d08SBarry Smith SNES_FAS *fas; 892ab8d36c9SPeter Brune PetscErrorCode ierr; 893ab8d36c9SPeter Brune SNES levelsnes; 89422d28d08SBarry Smith 895ab8d36c9SPeter Brune PetscFunctionBegin; 896ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 897ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 898ab8d36c9SPeter Brune *mat = fas->interpolate; 899ab8d36c9SPeter Brune PetscFunctionReturn(0); 900ab8d36c9SPeter Brune } 901ab8d36c9SPeter Brune 902ab8d36c9SPeter Brune #undef __FUNCT__ 903ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 904ab8d36c9SPeter Brune /*@ 905ab8d36c9SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 906ab8d36c9SPeter Brune from level l to l-1. 907ab8d36c9SPeter Brune 908ab8d36c9SPeter Brune Input Parameters: 909ab8d36c9SPeter Brune + snes - the multigrid context 910ab8d36c9SPeter Brune . mat - the restriction matrix 911ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 912ab8d36c9SPeter Brune 913ab8d36c9SPeter Brune Level: advanced 914ab8d36c9SPeter Brune 915ab8d36c9SPeter Brune Notes: 916ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation 917ab8d36c9SPeter Brune for the same level. 918ab8d36c9SPeter Brune 919ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 920ab8d36c9SPeter Brune out from the matrix size which one it is. 921ab8d36c9SPeter Brune 922ab8d36c9SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 923ab8d36c9SPeter Brune is used. 924ab8d36c9SPeter Brune 925ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 926ab8d36c9SPeter Brune 927ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 928ab8d36c9SPeter Brune @*/ 92922d28d08SBarry Smith PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) 93022d28d08SBarry Smith { 93122d28d08SBarry Smith SNES_FAS *fas; 932ab8d36c9SPeter Brune PetscErrorCode ierr; 933ab8d36c9SPeter Brune SNES levelsnes; 93422d28d08SBarry Smith 935ab8d36c9SPeter Brune PetscFunctionBegin; 936ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 937ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 938ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 939ab8d36c9SPeter Brune ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 9401aa26658SKarl Rupp 941ab8d36c9SPeter Brune fas->restrct = mat; 942ab8d36c9SPeter Brune PetscFunctionReturn(0); 943ab8d36c9SPeter Brune } 944ab8d36c9SPeter Brune 945ab8d36c9SPeter Brune #undef __FUNCT__ 946ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetRestriction" 947ab8d36c9SPeter Brune /*@ 948ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the 949ab8d36c9SPeter Brune restriction from l to the l-1th level 950ab8d36c9SPeter Brune 951ab8d36c9SPeter Brune Input Parameters: 952ab8d36c9SPeter Brune + snes - the multigrid context 953ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 954ab8d36c9SPeter Brune 955ab8d36c9SPeter Brune Output Parameters: 956ab8d36c9SPeter Brune . mat - the interpolation operator 957ab8d36c9SPeter Brune 958ab8d36c9SPeter Brune Level: advanced 959ab8d36c9SPeter Brune 960ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, restrict, level 961ab8d36c9SPeter Brune 962ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale() 963ab8d36c9SPeter Brune @*/ 96422d28d08SBarry Smith PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) 96522d28d08SBarry Smith { 96622d28d08SBarry Smith SNES_FAS *fas; 967ab8d36c9SPeter Brune PetscErrorCode ierr; 968ab8d36c9SPeter Brune SNES levelsnes; 96922d28d08SBarry Smith 970ab8d36c9SPeter Brune PetscFunctionBegin; 971ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 972ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 973ab8d36c9SPeter Brune *mat = fas->restrct; 974ab8d36c9SPeter Brune PetscFunctionReturn(0); 975ab8d36c9SPeter Brune } 976ab8d36c9SPeter Brune 977ab8d36c9SPeter Brune 978ab8d36c9SPeter Brune #undef __FUNCT__ 979ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInjection" 980ab8d36c9SPeter Brune /*@ 981ab8d36c9SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 982ab8d36c9SPeter Brune from level l to l-1. 983ab8d36c9SPeter Brune 984ab8d36c9SPeter Brune Input Parameters: 985ab8d36c9SPeter Brune + snes - the multigrid context 986ab8d36c9SPeter Brune . mat - the restriction matrix 987ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 988ab8d36c9SPeter Brune 989ab8d36c9SPeter Brune Level: advanced 990ab8d36c9SPeter Brune 991ab8d36c9SPeter Brune Notes: 992ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to 993ab8d36c9SPeter Brune project the solution instead. 994ab8d36c9SPeter Brune 995ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 996ab8d36c9SPeter Brune 997ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 998ab8d36c9SPeter Brune @*/ 99922d28d08SBarry Smith PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) 100022d28d08SBarry Smith { 100122d28d08SBarry Smith SNES_FAS *fas; 1002ab8d36c9SPeter Brune PetscErrorCode ierr; 1003ab8d36c9SPeter Brune SNES levelsnes; 100422d28d08SBarry Smith 1005ab8d36c9SPeter Brune PetscFunctionBegin; 1006ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1007ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1008ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1009ab8d36c9SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 10101aa26658SKarl Rupp 1011ab8d36c9SPeter Brune fas->inject = mat; 1012ab8d36c9SPeter Brune PetscFunctionReturn(0); 1013ab8d36c9SPeter Brune } 1014ab8d36c9SPeter Brune 1015ab8d36c9SPeter Brune 1016ab8d36c9SPeter Brune #undef __FUNCT__ 1017ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInjection" 1018ab8d36c9SPeter Brune /*@ 1019ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the 1020ab8d36c9SPeter Brune injection from l-1 to the lth level 1021ab8d36c9SPeter Brune 1022ab8d36c9SPeter Brune Input Parameters: 1023ab8d36c9SPeter Brune + snes - the multigrid context 1024ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 1025ab8d36c9SPeter Brune 1026ab8d36c9SPeter Brune Output Parameters: 1027ab8d36c9SPeter Brune . mat - the injection operator 1028ab8d36c9SPeter Brune 1029ab8d36c9SPeter Brune Level: advanced 1030ab8d36c9SPeter Brune 1031ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, restrict, level 1032ab8d36c9SPeter Brune 1033ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale() 1034ab8d36c9SPeter Brune @*/ 103522d28d08SBarry Smith PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) 103622d28d08SBarry Smith { 103722d28d08SBarry Smith SNES_FAS *fas; 1038ab8d36c9SPeter Brune PetscErrorCode ierr; 1039ab8d36c9SPeter Brune SNES levelsnes; 104022d28d08SBarry Smith 1041ab8d36c9SPeter Brune PetscFunctionBegin; 1042ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1043ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1044ab8d36c9SPeter Brune *mat = fas->inject; 1045ab8d36c9SPeter Brune PetscFunctionReturn(0); 1046ab8d36c9SPeter Brune } 1047ab8d36c9SPeter Brune 1048ab8d36c9SPeter Brune #undef __FUNCT__ 1049ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 1050ab8d36c9SPeter Brune /*@ 1051ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 1052ab8d36c9SPeter Brune operator from level l to l-1. 1053ab8d36c9SPeter Brune 1054ab8d36c9SPeter Brune Input Parameters: 1055ab8d36c9SPeter Brune + snes - the multigrid context 1056ab8d36c9SPeter Brune . rscale - the restriction scaling 1057ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 1058ab8d36c9SPeter Brune 1059ab8d36c9SPeter Brune Level: advanced 1060ab8d36c9SPeter Brune 1061ab8d36c9SPeter Brune Notes: 1062ab8d36c9SPeter Brune This is only used in the case that the injection is not set. 1063ab8d36c9SPeter Brune 1064ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 1065ab8d36c9SPeter Brune 1066ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1067ab8d36c9SPeter Brune @*/ 106822d28d08SBarry Smith PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) 106922d28d08SBarry Smith { 1070ab8d36c9SPeter Brune SNES_FAS *fas; 1071ab8d36c9SPeter Brune PetscErrorCode ierr; 1072ab8d36c9SPeter Brune SNES levelsnes; 107322d28d08SBarry Smith 1074ab8d36c9SPeter Brune PetscFunctionBegin; 1075ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1076ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1077ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 1078ab8d36c9SPeter Brune ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 10791aa26658SKarl Rupp 1080ab8d36c9SPeter Brune fas->rscale = rscale; 1081ab8d36c9SPeter Brune PetscFunctionReturn(0); 1082ab8d36c9SPeter Brune } 1083ab8d36c9SPeter Brune 1084ab8d36c9SPeter Brune #undef __FUNCT__ 1085ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmoother" 1086ab8d36c9SPeter Brune /*@ 1087ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level. 1088ab8d36c9SPeter Brune 1089ab8d36c9SPeter Brune Input Parameters: 1090ab8d36c9SPeter Brune + snes - the multigrid context 1091ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1092ab8d36c9SPeter Brune 1093ab8d36c9SPeter Brune Output Parameters: 1094ab8d36c9SPeter Brune smooth - the smoother 1095ab8d36c9SPeter Brune 1096ab8d36c9SPeter Brune Level: advanced 1097ab8d36c9SPeter Brune 1098ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1099ab8d36c9SPeter Brune 1100ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1101ab8d36c9SPeter Brune @*/ 110222d28d08SBarry Smith PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) 110322d28d08SBarry Smith { 1104ab8d36c9SPeter Brune SNES_FAS *fas; 1105ab8d36c9SPeter Brune PetscErrorCode ierr; 1106ab8d36c9SPeter Brune SNES levelsnes; 110722d28d08SBarry Smith 1108ab8d36c9SPeter Brune PetscFunctionBegin; 1109ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1110ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1111ab8d36c9SPeter Brune if (!fas->smoothd) { 11123de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1113ab8d36c9SPeter Brune } 1114ab8d36c9SPeter Brune *smooth = fas->smoothd; 1115ab8d36c9SPeter Brune PetscFunctionReturn(0); 1116ab8d36c9SPeter Brune } 1117ab8d36c9SPeter Brune 1118ab8d36c9SPeter Brune #undef __FUNCT__ 1119ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherDown" 1120ab8d36c9SPeter Brune /*@ 1121ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level. 1122ab8d36c9SPeter Brune 1123ab8d36c9SPeter Brune Input Parameters: 1124ab8d36c9SPeter Brune + snes - the multigrid context 1125ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1126ab8d36c9SPeter Brune 1127ab8d36c9SPeter Brune Output Parameters: 1128ab8d36c9SPeter Brune smooth - the smoother 1129ab8d36c9SPeter Brune 1130ab8d36c9SPeter Brune Level: advanced 1131ab8d36c9SPeter Brune 1132ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1133ab8d36c9SPeter Brune 1134ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1135ab8d36c9SPeter Brune @*/ 113622d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) 113722d28d08SBarry Smith { 1138ab8d36c9SPeter Brune SNES_FAS *fas; 1139ab8d36c9SPeter Brune PetscErrorCode ierr; 1140ab8d36c9SPeter Brune SNES levelsnes; 114122d28d08SBarry Smith 1142ab8d36c9SPeter Brune PetscFunctionBegin; 1143ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1144ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1145ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1146ab8d36c9SPeter Brune if (!fas->smoothd) { 11473de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1148ab8d36c9SPeter Brune } 1149ab8d36c9SPeter Brune if (!fas->smoothu) { 11503de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr); 1151ab8d36c9SPeter Brune } 1152ab8d36c9SPeter Brune *smooth = fas->smoothd; 1153ab8d36c9SPeter Brune PetscFunctionReturn(0); 1154ab8d36c9SPeter Brune } 1155ab8d36c9SPeter Brune 1156ab8d36c9SPeter Brune #undef __FUNCT__ 1157ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherUp" 1158ab8d36c9SPeter Brune /*@ 1159ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level. 1160ab8d36c9SPeter Brune 1161ab8d36c9SPeter Brune Input Parameters: 1162ab8d36c9SPeter Brune + snes - the multigrid context 1163ab8d36c9SPeter Brune - level - the level (0 is coarsest) 1164ab8d36c9SPeter Brune 1165ab8d36c9SPeter Brune Output Parameters: 1166ab8d36c9SPeter Brune smooth - the smoother 1167ab8d36c9SPeter Brune 1168ab8d36c9SPeter Brune Level: advanced 1169ab8d36c9SPeter Brune 1170ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1171ab8d36c9SPeter Brune 1172ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1173ab8d36c9SPeter Brune @*/ 117422d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) 117522d28d08SBarry Smith { 1176ab8d36c9SPeter Brune SNES_FAS *fas; 1177ab8d36c9SPeter Brune PetscErrorCode ierr; 1178ab8d36c9SPeter Brune SNES levelsnes; 117922d28d08SBarry Smith 1180ab8d36c9SPeter Brune PetscFunctionBegin; 1181ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1182ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1183ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1184ab8d36c9SPeter Brune if (!fas->smoothd) { 11853de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1186ab8d36c9SPeter Brune } 1187ab8d36c9SPeter Brune if (!fas->smoothu) { 11883de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr); 1189ab8d36c9SPeter Brune } 1190ab8d36c9SPeter Brune *smooth = fas->smoothu; 1191ab8d36c9SPeter Brune PetscFunctionReturn(0); 1192ab8d36c9SPeter Brune } 1193d6ad1212SPeter Brune 1194d6ad1212SPeter Brune #undef __FUNCT__ 1195d6ad1212SPeter Brune #define __FUNCT__ "SNESFASGetCoarseSolve" 1196d6ad1212SPeter Brune /*@ 1197d6ad1212SPeter Brune SNESFASGetCoarseSolve - Gets the coarsest solver. 1198d6ad1212SPeter Brune 1199d6ad1212SPeter Brune Input Parameters: 1200a3a80b83SMatthew G. Knepley . snes - the multigrid context 1201d6ad1212SPeter Brune 1202d6ad1212SPeter Brune Output Parameters: 1203a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver 1204d6ad1212SPeter Brune 1205d6ad1212SPeter Brune Level: advanced 1206d6ad1212SPeter Brune 1207d6ad1212SPeter Brune .keywords: FAS, MG, get, multigrid, solver, coarse 1208d6ad1212SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1209d6ad1212SPeter Brune @*/ 1210a3a80b83SMatthew G. Knepley PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse) 121122d28d08SBarry Smith { 1212d6ad1212SPeter Brune SNES_FAS *fas; 1213d6ad1212SPeter Brune PetscErrorCode ierr; 1214d6ad1212SPeter Brune SNES levelsnes; 121522d28d08SBarry Smith 1216d6ad1212SPeter Brune PetscFunctionBegin; 1217d6ad1212SPeter Brune ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr); 1218d6ad1212SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1219d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1220d6ad1212SPeter Brune if (!fas->smoothd) { 12213de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1222d6ad1212SPeter Brune } 1223a3a80b83SMatthew G. Knepley *coarse = fas->smoothd; 1224d6ad1212SPeter Brune PetscFunctionReturn(0); 1225d6ad1212SPeter Brune } 1226928e959bSPeter Brune 1227928e959bSPeter Brune #undef __FUNCT__ 1228928e959bSPeter Brune #define __FUNCT__ "SNESFASFullSetDownSweep" 1229928e959bSPeter Brune /*@ 1230928e959bSPeter Brune SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS 1231928e959bSPeter Brune 1232928e959bSPeter Brune Logically Collective on SNES 1233928e959bSPeter Brune 1234928e959bSPeter Brune Input Parameters: 1235928e959bSPeter Brune + snes - the multigrid context 1236928e959bSPeter Brune - swp - whether to downsweep or not 1237928e959bSPeter Brune 1238928e959bSPeter Brune Options Database Key: 1239928e959bSPeter Brune . -snes_fas_full_downsweep - Sets number of pre-smoothing steps 1240928e959bSPeter Brune 1241928e959bSPeter Brune Level: advanced 1242928e959bSPeter Brune 1243928e959bSPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 1244928e959bSPeter Brune 1245928e959bSPeter Brune .seealso: SNESFASSetNumberSmoothUp() 1246928e959bSPeter Brune @*/ 1247928e959bSPeter Brune PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp) 1248928e959bSPeter Brune { 1249928e959bSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 1250928e959bSPeter Brune PetscErrorCode ierr = 0; 1251928e959bSPeter Brune 1252928e959bSPeter Brune PetscFunctionBegin; 1253928e959bSPeter Brune fas->full_downsweep = swp; 1254928e959bSPeter Brune if (fas->next) { 1255928e959bSPeter Brune ierr = SNESFASFullSetDownSweep(fas->next,swp);CHKERRQ(ierr); 1256928e959bSPeter Brune } 1257928e959bSPeter Brune PetscFunctionReturn(0); 1258928e959bSPeter Brune } 1259