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 77*2bf68e3eSBarry Smith problems on smaller sets of processors. 78ab8d36c9SPeter Brune 79ab8d36c9SPeter Brune Level: intermediate 80ab8d36c9SPeter Brune 81ab8d36c9SPeter Brune Notes: 82ab8d36c9SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 83ab8d36c9SPeter Brune for setting the level options rather than the -fas_coarse prefix. 84ab8d36c9SPeter Brune 85ab8d36c9SPeter Brune .keywords: FAS, MG, set, levels, multigrid 86ab8d36c9SPeter Brune 87ab8d36c9SPeter Brune .seealso: SNESFASGetLevels() 88ab8d36c9SPeter Brune @*/ 8922d28d08SBarry Smith PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) 9022d28d08SBarry Smith { 91ab8d36c9SPeter Brune PetscErrorCode ierr; 92ab8d36c9SPeter Brune PetscInt i; 93ab8d36c9SPeter Brune const char *optionsprefix; 94ab8d36c9SPeter Brune char tprefix[128]; 95ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 96ab8d36c9SPeter Brune SNES prevsnes; 97ab8d36c9SPeter Brune MPI_Comm comm; 9822d28d08SBarry Smith 99ab8d36c9SPeter Brune PetscFunctionBegin; 100ce94432eSBarry Smith ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr); 101ab8d36c9SPeter Brune if (levels == fas->levels) { 10222d28d08SBarry Smith if (!comms) PetscFunctionReturn(0); 103ab8d36c9SPeter Brune } 104ab8d36c9SPeter Brune /* user has changed the number of levels; reset */ 105ab8d36c9SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 106ab8d36c9SPeter Brune /* destroy any coarser levels if necessary */ 107ab8d36c9SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 1080298fd71SBarry Smith fas->next = NULL; 1090298fd71SBarry Smith fas->previous = NULL; 110ab8d36c9SPeter Brune prevsnes = snes; 111ab8d36c9SPeter Brune /* setup the finest level */ 112ab8d36c9SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 113ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1);CHKERRQ(ierr); 114ab8d36c9SPeter Brune for (i = levels - 1; i >= 0; i--) { 115ab8d36c9SPeter Brune if (comms) comm = comms[i]; 116ab8d36c9SPeter Brune fas->level = i; 117ab8d36c9SPeter Brune fas->levels = levels; 118ab8d36c9SPeter Brune fas->fine = snes; 1190298fd71SBarry Smith fas->next = NULL; 120ab8d36c9SPeter Brune if (i > 0) { 121ab8d36c9SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 122e964f0dbSPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 123ab8d36c9SPeter Brune sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level); 124ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(fas->next,optionsprefix);CHKERRQ(ierr); 125e964f0dbSPeter Brune ierr = SNESAppendOptionsPrefix(fas->next,tprefix);CHKERRQ(ierr); 126ab8d36c9SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 127ab8d36c9SPeter Brune ierr = SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);CHKERRQ(ierr); 128ab8d36c9SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 129ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1);CHKERRQ(ierr); 1301aa26658SKarl Rupp 131ab8d36c9SPeter Brune ((SNES_FAS*)fas->next->data)->previous = prevsnes; 1321aa26658SKarl Rupp 133ab8d36c9SPeter Brune prevsnes = fas->next; 134ab8d36c9SPeter Brune fas = (SNES_FAS*)prevsnes->data; 135ab8d36c9SPeter Brune } 136ab8d36c9SPeter Brune } 137ab8d36c9SPeter Brune PetscFunctionReturn(0); 138ab8d36c9SPeter Brune } 139ab8d36c9SPeter Brune 140ab8d36c9SPeter Brune 141ab8d36c9SPeter Brune #undef __FUNCT__ 142ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 143ab8d36c9SPeter Brune /*@ 144ab8d36c9SPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 145ab8d36c9SPeter Brune 146ab8d36c9SPeter Brune Input Parameter: 147ab8d36c9SPeter Brune . snes - the nonlinear solver context 148ab8d36c9SPeter Brune 149ab8d36c9SPeter Brune Output parameter: 150ab8d36c9SPeter Brune . levels - the number of levels 151ab8d36c9SPeter Brune 152ab8d36c9SPeter Brune Level: advanced 153ab8d36c9SPeter Brune 154ab8d36c9SPeter Brune .keywords: MG, get, levels, multigrid 155ab8d36c9SPeter Brune 156ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 157ab8d36c9SPeter Brune @*/ 15822d28d08SBarry Smith PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels) 15922d28d08SBarry Smith { 160ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 1615fd66863SKarl Rupp 162ab8d36c9SPeter Brune PetscFunctionBegin; 163ab8d36c9SPeter Brune *levels = fas->levels; 164ab8d36c9SPeter Brune PetscFunctionReturn(0); 165ab8d36c9SPeter Brune } 166ab8d36c9SPeter Brune 167ab8d36c9SPeter Brune 168ab8d36c9SPeter Brune #undef __FUNCT__ 169ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetCycleSNES" 170ab8d36c9SPeter Brune /*@ 171ab8d36c9SPeter Brune SNESFASGetCycleSNES - Gets the SNES corresponding to a particular 172ab8d36c9SPeter Brune level of the FAS hierarchy. 173ab8d36c9SPeter Brune 174ab8d36c9SPeter Brune Input Parameters: 175ab8d36c9SPeter Brune + snes - the multigrid context 176ab8d36c9SPeter Brune level - the level to get 177ab8d36c9SPeter Brune - lsnes - whether to use the nonlinear smoother or not 178ab8d36c9SPeter Brune 179ab8d36c9SPeter Brune Level: advanced 180ab8d36c9SPeter Brune 181ab8d36c9SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 182ab8d36c9SPeter Brune 183ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 184ab8d36c9SPeter Brune @*/ 18522d28d08SBarry Smith PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes) 18622d28d08SBarry Smith { 187ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 188ab8d36c9SPeter Brune PetscInt i; 189ab8d36c9SPeter Brune 190ab8d36c9SPeter Brune PetscFunctionBegin; 191ce94432eSBarry 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); 192ce94432eSBarry 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); 193ab8d36c9SPeter Brune 194ab8d36c9SPeter Brune *lsnes = snes; 195ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) { 196ab8d36c9SPeter Brune *lsnes = fas->next; 197ab8d36c9SPeter Brune fas = (SNES_FAS*)(*lsnes)->data; 198ab8d36c9SPeter Brune } 199ce94432eSBarry Smith if (fas->level != level) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 200ab8d36c9SPeter Brune PetscFunctionReturn(0); 201ab8d36c9SPeter Brune } 202ab8d36c9SPeter Brune 203ab8d36c9SPeter Brune #undef __FUNCT__ 204ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 205ab8d36c9SPeter Brune /*@ 206ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 207ab8d36c9SPeter Brune use on all levels. 208ab8d36c9SPeter Brune 209ab8d36c9SPeter Brune Logically Collective on SNES 210ab8d36c9SPeter Brune 211ab8d36c9SPeter Brune Input Parameters: 212ab8d36c9SPeter Brune + snes - the multigrid context 213ab8d36c9SPeter Brune - n - the number of smoothing steps 214ab8d36c9SPeter Brune 215ab8d36c9SPeter Brune Options Database Key: 216ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 217ab8d36c9SPeter Brune 218ab8d36c9SPeter Brune Level: advanced 219ab8d36c9SPeter Brune 220ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 221ab8d36c9SPeter Brune 222ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 223ab8d36c9SPeter Brune @*/ 22422d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) 22522d28d08SBarry Smith { 226ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 22722d28d08SBarry Smith PetscErrorCode ierr; 22822d28d08SBarry Smith 229ab8d36c9SPeter Brune PetscFunctionBegin; 230ab8d36c9SPeter Brune fas->max_up_it = n; 231656ede7eSPeter Brune if (!fas->smoothu && fas->level != 0) { 23222d28d08SBarry Smith ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 233ab8d36c9SPeter Brune } 23422d28d08SBarry Smith if (fas->smoothu) { 23522d28d08SBarry Smith ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr); 23622d28d08SBarry Smith } 237ab8d36c9SPeter Brune if (fas->next) { 238ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 239ab8d36c9SPeter Brune } 240ab8d36c9SPeter Brune PetscFunctionReturn(0); 241ab8d36c9SPeter Brune } 242ab8d36c9SPeter Brune 243ab8d36c9SPeter Brune #undef __FUNCT__ 244ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 245ab8d36c9SPeter Brune /*@ 246ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 247ab8d36c9SPeter Brune use on all levels. 248ab8d36c9SPeter Brune 249ab8d36c9SPeter Brune Logically Collective on SNES 250ab8d36c9SPeter Brune 251ab8d36c9SPeter Brune Input Parameters: 252ab8d36c9SPeter Brune + snes - the multigrid context 253ab8d36c9SPeter Brune - n - the number of smoothing steps 254ab8d36c9SPeter Brune 255ab8d36c9SPeter Brune Options Database Key: 256ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 257ab8d36c9SPeter Brune 258ab8d36c9SPeter Brune Level: advanced 259ab8d36c9SPeter Brune 260ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 261ab8d36c9SPeter Brune 262ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 263ab8d36c9SPeter Brune @*/ 26422d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) 26522d28d08SBarry Smith { 266ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 267ab8d36c9SPeter Brune PetscErrorCode ierr = 0; 26822d28d08SBarry Smith 269ab8d36c9SPeter Brune PetscFunctionBegin; 270ab8d36c9SPeter Brune if (!fas->smoothd) { 27122d28d08SBarry Smith ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 272ab8d36c9SPeter Brune } 273ab8d36c9SPeter Brune ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr); 2741aa26658SKarl Rupp 275ab8d36c9SPeter Brune fas->max_down_it = n; 276ab8d36c9SPeter Brune if (fas->next) { 277ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 278ab8d36c9SPeter Brune } 279ab8d36c9SPeter Brune PetscFunctionReturn(0); 280ab8d36c9SPeter Brune } 281ab8d36c9SPeter Brune 282ab8d36c9SPeter Brune 283ab8d36c9SPeter Brune #undef __FUNCT__ 28487f44e3fSPeter Brune #define __FUNCT__ "SNESFASSetContinuation" 28587f44e3fSPeter Brune /*@ 28687f44e3fSPeter Brune SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep 28787f44e3fSPeter Brune 28887f44e3fSPeter Brune Logically Collective on SNES 28987f44e3fSPeter Brune 29087f44e3fSPeter Brune Input Parameters: 29187f44e3fSPeter Brune + snes - the multigrid context 29287f44e3fSPeter Brune - n - the number of smoothing steps 29387f44e3fSPeter Brune 29487f44e3fSPeter Brune Options Database Key: 29587f44e3fSPeter Brune . -snes_fas_continuation - sets continuation to true 29687f44e3fSPeter Brune 29787f44e3fSPeter Brune Level: advanced 29887f44e3fSPeter Brune 29987f44e3fSPeter Brune Notes: This sets the prefix on the upsweep smoothers to -fas_continuation 30087f44e3fSPeter Brune 30187f44e3fSPeter Brune .keywords: FAS, MG, smoother, continuation 30287f44e3fSPeter Brune 30387f44e3fSPeter Brune .seealso: SNESFAS 30487f44e3fSPeter Brune @*/ 30587f44e3fSPeter Brune PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation) 30687f44e3fSPeter Brune { 30787f44e3fSPeter Brune const char *optionsprefix; 30887f44e3fSPeter Brune char tprefix[128]; 30987f44e3fSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 31087f44e3fSPeter Brune PetscErrorCode ierr = 0; 31187f44e3fSPeter Brune 31287f44e3fSPeter Brune PetscFunctionBegin; 31387f44e3fSPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 31487f44e3fSPeter Brune if (!fas->smoothu) { 31587f44e3fSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 31687f44e3fSPeter Brune } 31787f44e3fSPeter Brune sprintf(tprefix,"fas_levels_continuation_"); 31887f44e3fSPeter Brune ierr = SNESSetOptionsPrefix(fas->smoothu, optionsprefix);CHKERRQ(ierr); 31987f44e3fSPeter Brune ierr = SNESAppendOptionsPrefix(fas->smoothu, tprefix);CHKERRQ(ierr); 32087f44e3fSPeter Brune ierr = SNESSetType(fas->smoothu,SNESNEWTONLS);CHKERRQ(ierr); 32187f44e3fSPeter Brune ierr = SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100);CHKERRQ(ierr); 32287f44e3fSPeter Brune fas->continuation = continuation; 32387f44e3fSPeter Brune if (fas->next) { 32487f44e3fSPeter Brune ierr = SNESFASSetContinuation(fas->next,continuation);CHKERRQ(ierr); 32587f44e3fSPeter Brune } 32687f44e3fSPeter Brune PetscFunctionReturn(0); 32787f44e3fSPeter Brune } 32887f44e3fSPeter Brune 32987f44e3fSPeter Brune 33087f44e3fSPeter Brune #undef __FUNCT__ 331ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetCycles" 332ab8d36c9SPeter Brune /*@ 333ab8d36c9SPeter Brune SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 334ab8d36c9SPeter Brune complicated cycling. 335ab8d36c9SPeter Brune 336ab8d36c9SPeter Brune Logically Collective on SNES 337ab8d36c9SPeter Brune 338ab8d36c9SPeter Brune Input Parameters: 339ab8d36c9SPeter Brune + snes - the multigrid context 340ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 341ab8d36c9SPeter Brune 342ab8d36c9SPeter Brune Options Database Key: 343e1bc860dSBarry Smith . -snes_fas_cycles 1 or 2 344ab8d36c9SPeter Brune 345ab8d36c9SPeter Brune Level: advanced 346ab8d36c9SPeter Brune 347ab8d36c9SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 348ab8d36c9SPeter Brune 349ab8d36c9SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 350ab8d36c9SPeter Brune @*/ 35122d28d08SBarry Smith PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) 35222d28d08SBarry Smith { 353ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 354ab8d36c9SPeter Brune PetscErrorCode ierr; 355ab8d36c9SPeter Brune PetscBool isFine; 35622d28d08SBarry Smith 357ab8d36c9SPeter Brune PetscFunctionBegin; 35822d28d08SBarry Smith ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 3591aa26658SKarl Rupp 360ab8d36c9SPeter Brune fas->n_cycles = cycles; 3611aa26658SKarl Rupp if (!isFine) { 362ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 3631aa26658SKarl Rupp } 364ab8d36c9SPeter Brune if (fas->next) { 365ab8d36c9SPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 366ab8d36c9SPeter Brune } 367ab8d36c9SPeter Brune PetscFunctionReturn(0); 368ab8d36c9SPeter Brune } 369ab8d36c9SPeter Brune 370c8c899caSPeter Brune 371c8c899caSPeter Brune #undef __FUNCT__ 372c8c899caSPeter Brune #define __FUNCT__ "SNESFASSetMonitor" 373c8c899caSPeter Brune /*@ 374c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring 375c8c899caSPeter Brune 376c8c899caSPeter Brune Logically Collective on SNES 377c8c899caSPeter Brune 378c8c899caSPeter Brune Input Parameters: 379c8c899caSPeter Brune + snes - the FAS context 380d142ab34SLawrence Mitchell . vf - viewer and format structure (may be NULL if flg is FALSE) 381c8c899caSPeter Brune - flg - monitor or not 382c8c899caSPeter Brune 383c8c899caSPeter Brune Level: advanced 384c8c899caSPeter Brune 385c8c899caSPeter Brune .keywords: FAS, monitor 386c8c899caSPeter Brune 387c8c899caSPeter Brune .seealso: SNESFASSetCyclesOnLevel() 388c8c899caSPeter Brune @*/ 389d142ab34SLawrence Mitchell PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg) 39022d28d08SBarry Smith { 391c8c899caSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 392c8c899caSPeter Brune PetscErrorCode ierr; 393c8c899caSPeter Brune PetscBool isFine; 394c8c899caSPeter Brune PetscInt i, levels = fas->levels; 395c8c899caSPeter Brune SNES levelsnes; 39622d28d08SBarry Smith 397c8c899caSPeter Brune PetscFunctionBegin; 398c8c899caSPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 399c8c899caSPeter Brune if (isFine) { 400c8c899caSPeter Brune for (i = 0; i < levels; i++) { 40122d28d08SBarry Smith ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 402c8c899caSPeter Brune fas = (SNES_FAS*)levelsnes->data; 403c8c899caSPeter Brune if (flg) { 404c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */ 405c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 406d142ab34SLawrence Mitchell /* Only register destroy on finest level */ 407d142ab34SLawrence Mitchell ierr = SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL));CHKERRQ(ierr); 4081aa26658SKarl Rupp } else if (i != fas->levels - 1) { 409c8c899caSPeter Brune /* unset the monitors on the coarse levels */ 410c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 411c8c899caSPeter Brune } 412c8c899caSPeter Brune } 413c8c899caSPeter Brune } 414c8c899caSPeter Brune PetscFunctionReturn(0); 415c8c899caSPeter Brune } 416c8c899caSPeter Brune 417ab8d36c9SPeter Brune #undef __FUNCT__ 4180dd27c6cSPeter Brune #define __FUNCT__ "SNESFASSetLog" 4190dd27c6cSPeter Brune /*@ 4200dd27c6cSPeter Brune SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels 4210dd27c6cSPeter Brune 4220dd27c6cSPeter Brune Logically Collective on SNES 4230dd27c6cSPeter Brune 4240dd27c6cSPeter Brune Input Parameters: 4250dd27c6cSPeter Brune + snes - the FAS context 4260dd27c6cSPeter Brune - flg - monitor or not 4270dd27c6cSPeter Brune 4280dd27c6cSPeter Brune Level: advanced 4290dd27c6cSPeter Brune 4300dd27c6cSPeter Brune .keywords: FAS, logging 4310dd27c6cSPeter Brune 4320dd27c6cSPeter Brune .seealso: SNESFASSetMonitor() 4330dd27c6cSPeter Brune @*/ 4340dd27c6cSPeter Brune PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) 4350dd27c6cSPeter Brune { 4360dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 4370dd27c6cSPeter Brune PetscErrorCode ierr; 4380dd27c6cSPeter Brune PetscBool isFine; 4390dd27c6cSPeter Brune PetscInt i, levels = fas->levels; 4400dd27c6cSPeter Brune SNES levelsnes; 4410dd27c6cSPeter Brune char eventname[128]; 4420dd27c6cSPeter Brune 4430dd27c6cSPeter Brune PetscFunctionBegin; 4440dd27c6cSPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 4450dd27c6cSPeter Brune if (isFine) { 4460dd27c6cSPeter Brune for (i = 0; i < levels; i++) { 4470dd27c6cSPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 4480dd27c6cSPeter Brune fas = (SNES_FAS*)levelsnes->data; 4490dd27c6cSPeter Brune if (flg) { 4500dd27c6cSPeter Brune sprintf(eventname,"FASSetup %d",(int)i); 4510dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr); 4520dd27c6cSPeter Brune sprintf(eventname,"FASSmooth %d",(int)i); 4530dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr); 4540dd27c6cSPeter Brune sprintf(eventname,"FASResid %d",(int)i); 4550dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr); 4560dd27c6cSPeter Brune if (i) { 4570dd27c6cSPeter Brune sprintf(eventname,"FASInterp %d",(int)i); 4580dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr); 4590dd27c6cSPeter Brune } 4600dd27c6cSPeter Brune } else { 4610298fd71SBarry Smith fas->eventsmoothsetup = 0; 4620298fd71SBarry Smith fas->eventsmoothsolve = 0; 4630298fd71SBarry Smith fas->eventresidual = 0; 4640298fd71SBarry Smith fas->eventinterprestrict = 0; 4650dd27c6cSPeter Brune } 4660dd27c6cSPeter Brune } 4670dd27c6cSPeter Brune } 4680dd27c6cSPeter Brune PetscFunctionReturn(0); 4690dd27c6cSPeter Brune } 4700dd27c6cSPeter Brune 4710dd27c6cSPeter Brune #undef __FUNCT__ 472ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleCreateSmoother_Private" 473ab8d36c9SPeter Brune /* 474ab8d36c9SPeter Brune Creates the default smoother type. 475ab8d36c9SPeter Brune 47604d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level. 477ab8d36c9SPeter Brune 478ab8d36c9SPeter Brune */ 47922d28d08SBarry Smith PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) 48022d28d08SBarry Smith { 481ab8d36c9SPeter Brune SNES_FAS *fas; 482ab8d36c9SPeter Brune const char *optionsprefix; 483ab8d36c9SPeter Brune char tprefix[128]; 484ab8d36c9SPeter Brune PetscErrorCode ierr; 485ab8d36c9SPeter Brune SNES nsmooth; 48622d28d08SBarry Smith 487ab8d36c9SPeter Brune PetscFunctionBegin; 488ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 489ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 490ab8d36c9SPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 491ab8d36c9SPeter Brune /* create the default smoother */ 492ce94432eSBarry Smith ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);CHKERRQ(ierr); 493ab8d36c9SPeter Brune if (fas->level == 0) { 494ab8d36c9SPeter Brune sprintf(tprefix,"fas_coarse_"); 495ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 496ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 49704d7464bSBarry Smith ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr); 498e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr); 499ab8d36c9SPeter Brune } else { 500ab8d36c9SPeter Brune sprintf(tprefix,"fas_levels_%d_",(int)fas->level); 501ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 502ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 503ab8d36c9SPeter Brune ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr); 504e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr); 505ab8d36c9SPeter Brune } 506ab8d36c9SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 5073bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);CHKERRQ(ierr); 508f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr); 509ab83eea4SMatthew G. Knepley ierr = PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level);CHKERRQ(ierr); 510ab8d36c9SPeter Brune *smooth = nsmooth; 511ab8d36c9SPeter Brune PetscFunctionReturn(0); 512ab8d36c9SPeter Brune } 513ab8d36c9SPeter Brune 514ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */ 515ab8d36c9SPeter Brune 516ab8d36c9SPeter Brune #undef __FUNCT__ 517ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleSetCycles" 518ab8d36c9SPeter Brune /*@ 519ab8d36c9SPeter Brune SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 520ab8d36c9SPeter Brune 521ab8d36c9SPeter Brune Logically Collective on SNES 522ab8d36c9SPeter Brune 523ab8d36c9SPeter Brune Input Parameters: 524ab8d36c9SPeter Brune + snes - the multigrid context 525ab8d36c9SPeter Brune . level - the level to set the number of cycles on 526ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 527ab8d36c9SPeter Brune 528ab8d36c9SPeter Brune Level: advanced 529ab8d36c9SPeter Brune 530ab8d36c9SPeter Brune .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid 531ab8d36c9SPeter Brune 532ab8d36c9SPeter Brune .seealso: SNESFASSetCycles() 533ab8d36c9SPeter Brune @*/ 53422d28d08SBarry Smith PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) 53522d28d08SBarry Smith { 536ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 537ab8d36c9SPeter Brune PetscErrorCode ierr; 53822d28d08SBarry Smith 539ab8d36c9SPeter Brune PetscFunctionBegin; 540ab8d36c9SPeter Brune fas->n_cycles = cycles; 541ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 542ab8d36c9SPeter Brune PetscFunctionReturn(0); 543ab8d36c9SPeter Brune } 544ab8d36c9SPeter Brune 545ab8d36c9SPeter Brune 546ab8d36c9SPeter Brune #undef __FUNCT__ 547ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmoother" 548ab8d36c9SPeter Brune /*@ 549ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 550ab8d36c9SPeter Brune 551ab8d36c9SPeter Brune Logically Collective on SNES 552ab8d36c9SPeter Brune 553ab8d36c9SPeter Brune Input Parameters: 554ab8d36c9SPeter Brune . snes - the multigrid context 555ab8d36c9SPeter Brune 556ab8d36c9SPeter Brune Output Parameters: 557ab8d36c9SPeter Brune . smooth - the smoother 558ab8d36c9SPeter Brune 559ab8d36c9SPeter Brune Level: advanced 560ab8d36c9SPeter Brune 561ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 562ab8d36c9SPeter Brune 563ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown() 564ab8d36c9SPeter Brune @*/ 565ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 566ab8d36c9SPeter Brune { 567ab8d36c9SPeter Brune SNES_FAS *fas; 5685fd66863SKarl Rupp 569ab8d36c9SPeter Brune PetscFunctionBegin; 570ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 571ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 572ab8d36c9SPeter Brune *smooth = fas->smoothd; 573ab8d36c9SPeter Brune PetscFunctionReturn(0); 574ab8d36c9SPeter Brune } 575ab8d36c9SPeter Brune #undef __FUNCT__ 576ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherUp" 577ab8d36c9SPeter Brune /*@ 578ab8d36c9SPeter Brune SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 579ab8d36c9SPeter Brune 580ab8d36c9SPeter Brune Logically Collective on SNES 581ab8d36c9SPeter Brune 582ab8d36c9SPeter Brune Input Parameters: 583ab8d36c9SPeter Brune . snes - the multigrid context 584ab8d36c9SPeter Brune 585ab8d36c9SPeter Brune Output Parameters: 586ab8d36c9SPeter Brune . smoothu - the smoother 587ab8d36c9SPeter Brune 588ab8d36c9SPeter Brune Notes: 589ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent 590ab8d36c9SPeter Brune default behavior in the process of the solve. 591ab8d36c9SPeter Brune 592ab8d36c9SPeter Brune Level: advanced 593ab8d36c9SPeter Brune 594ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 595ab8d36c9SPeter Brune 596ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown() 597ab8d36c9SPeter Brune @*/ 598ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 599ab8d36c9SPeter Brune { 600ab8d36c9SPeter Brune SNES_FAS *fas; 6015fd66863SKarl Rupp 602ab8d36c9SPeter Brune PetscFunctionBegin; 603ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 604ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 6051aa26658SKarl Rupp if (!fas->smoothu) *smoothu = fas->smoothd; 6061aa26658SKarl Rupp else *smoothu = fas->smoothu; 607ab8d36c9SPeter Brune PetscFunctionReturn(0); 608ab8d36c9SPeter Brune } 609ab8d36c9SPeter Brune 610ab8d36c9SPeter Brune #undef __FUNCT__ 611ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherDown" 612ab8d36c9SPeter Brune /*@ 613ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 614ab8d36c9SPeter Brune 615ab8d36c9SPeter Brune Logically Collective on SNES 616ab8d36c9SPeter Brune 617ab8d36c9SPeter Brune Input Parameters: 618ab8d36c9SPeter Brune . snes - the multigrid context 619ab8d36c9SPeter Brune 620ab8d36c9SPeter Brune Output Parameters: 621ab8d36c9SPeter Brune . smoothd - the smoother 622ab8d36c9SPeter Brune 623ab8d36c9SPeter Brune Level: advanced 624ab8d36c9SPeter Brune 625ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 626ab8d36c9SPeter Brune 627ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 628ab8d36c9SPeter Brune @*/ 629ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 630ab8d36c9SPeter Brune { 631ab8d36c9SPeter Brune SNES_FAS *fas; 6325fd66863SKarl Rupp 633ab8d36c9SPeter Brune PetscFunctionBegin; 634ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 635ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 636ab8d36c9SPeter Brune *smoothd = fas->smoothd; 637ab8d36c9SPeter Brune PetscFunctionReturn(0); 638ab8d36c9SPeter Brune } 639ab8d36c9SPeter Brune 640ab8d36c9SPeter Brune 641ab8d36c9SPeter Brune #undef __FUNCT__ 642ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetCorrection" 643ab8d36c9SPeter Brune /*@ 644ab8d36c9SPeter Brune SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 645ab8d36c9SPeter Brune 646ab8d36c9SPeter Brune Logically Collective on SNES 647ab8d36c9SPeter Brune 648ab8d36c9SPeter Brune Input Parameters: 649ab8d36c9SPeter Brune . snes - the multigrid context 650ab8d36c9SPeter Brune 651ab8d36c9SPeter Brune Output Parameters: 652ab8d36c9SPeter Brune . correction - the coarse correction on this level 653ab8d36c9SPeter Brune 654ab8d36c9SPeter Brune Notes: 6550298fd71SBarry Smith Returns NULL on the coarsest level. 656ab8d36c9SPeter Brune 657ab8d36c9SPeter Brune Level: advanced 658ab8d36c9SPeter Brune 659ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 660ab8d36c9SPeter Brune 661ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 662ab8d36c9SPeter Brune @*/ 663ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 664ab8d36c9SPeter Brune { 665ab8d36c9SPeter Brune SNES_FAS *fas; 6665fd66863SKarl Rupp 667ab8d36c9SPeter Brune PetscFunctionBegin; 668ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 669ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 670ab8d36c9SPeter Brune *correction = fas->next; 671ab8d36c9SPeter Brune PetscFunctionReturn(0); 672ab8d36c9SPeter Brune } 673ab8d36c9SPeter Brune 674ab8d36c9SPeter Brune #undef __FUNCT__ 675ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInterpolation" 676ab8d36c9SPeter Brune /*@ 677ab8d36c9SPeter Brune SNESFASCycleGetInterpolation - Gets the interpolation on this level 678ab8d36c9SPeter Brune 679ab8d36c9SPeter Brune Logically Collective on SNES 680ab8d36c9SPeter Brune 681ab8d36c9SPeter Brune Input Parameters: 682ab8d36c9SPeter Brune . snes - the multigrid context 683ab8d36c9SPeter Brune 684ab8d36c9SPeter Brune Output Parameters: 685ab8d36c9SPeter Brune . mat - the interpolation operator on this level 686ab8d36c9SPeter Brune 687ab8d36c9SPeter Brune Level: developer 688ab8d36c9SPeter Brune 689ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 690ab8d36c9SPeter Brune 691ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 692ab8d36c9SPeter Brune @*/ 693ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 694ab8d36c9SPeter Brune { 695ab8d36c9SPeter Brune SNES_FAS *fas; 6965fd66863SKarl Rupp 697ab8d36c9SPeter Brune PetscFunctionBegin; 698ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 699ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 700ab8d36c9SPeter Brune *mat = fas->interpolate; 701ab8d36c9SPeter Brune PetscFunctionReturn(0); 702ab8d36c9SPeter Brune } 703ab8d36c9SPeter Brune 704ab8d36c9SPeter Brune 705ab8d36c9SPeter Brune #undef __FUNCT__ 706ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRestriction" 707ab8d36c9SPeter Brune /*@ 708ab8d36c9SPeter Brune SNESFASCycleGetRestriction - Gets the restriction on this level 709ab8d36c9SPeter Brune 710ab8d36c9SPeter Brune Logically Collective on SNES 711ab8d36c9SPeter Brune 712ab8d36c9SPeter Brune Input Parameters: 713ab8d36c9SPeter Brune . snes - the multigrid context 714ab8d36c9SPeter Brune 715ab8d36c9SPeter Brune Output Parameters: 716ab8d36c9SPeter Brune . mat - the restriction operator on this level 717ab8d36c9SPeter Brune 718ab8d36c9SPeter Brune Level: developer 719ab8d36c9SPeter Brune 720ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 721ab8d36c9SPeter Brune 722ab8d36c9SPeter Brune .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation() 723ab8d36c9SPeter Brune @*/ 724ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 725ab8d36c9SPeter Brune { 726ab8d36c9SPeter Brune SNES_FAS *fas; 7275fd66863SKarl Rupp 728ab8d36c9SPeter Brune PetscFunctionBegin; 729ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 730ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 731ab8d36c9SPeter Brune *mat = fas->restrct; 732ab8d36c9SPeter Brune PetscFunctionReturn(0); 733ab8d36c9SPeter Brune } 734ab8d36c9SPeter Brune 735ab8d36c9SPeter Brune 736ab8d36c9SPeter Brune #undef __FUNCT__ 737ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInjection" 738ab8d36c9SPeter Brune /*@ 739ab8d36c9SPeter Brune SNESFASCycleGetInjection - Gets the injection on this level 740ab8d36c9SPeter Brune 741ab8d36c9SPeter Brune Logically Collective on SNES 742ab8d36c9SPeter Brune 743ab8d36c9SPeter Brune Input Parameters: 744ab8d36c9SPeter Brune . snes - the multigrid context 745ab8d36c9SPeter Brune 746ab8d36c9SPeter Brune Output Parameters: 747ab8d36c9SPeter Brune . mat - the restriction operator on this level 748ab8d36c9SPeter Brune 749ab8d36c9SPeter Brune Level: developer 750ab8d36c9SPeter Brune 751ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 752ab8d36c9SPeter Brune 753ab8d36c9SPeter Brune .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction() 754ab8d36c9SPeter Brune @*/ 755ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 756ab8d36c9SPeter Brune { 757ab8d36c9SPeter Brune SNES_FAS *fas; 7585fd66863SKarl Rupp 759ab8d36c9SPeter Brune PetscFunctionBegin; 760ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 761ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 762ab8d36c9SPeter Brune *mat = fas->inject; 763ab8d36c9SPeter Brune PetscFunctionReturn(0); 764ab8d36c9SPeter Brune } 765ab8d36c9SPeter Brune 766ab8d36c9SPeter Brune #undef __FUNCT__ 767ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRScale" 768ab8d36c9SPeter Brune /*@ 769ab8d36c9SPeter Brune SNESFASCycleGetRScale - Gets the injection on this level 770ab8d36c9SPeter Brune 771ab8d36c9SPeter Brune Logically Collective on SNES 772ab8d36c9SPeter Brune 773ab8d36c9SPeter Brune Input Parameters: 774ab8d36c9SPeter Brune . snes - the multigrid context 775ab8d36c9SPeter Brune 776ab8d36c9SPeter Brune Output Parameters: 777ab8d36c9SPeter Brune . mat - the restriction operator on this level 778ab8d36c9SPeter Brune 779ab8d36c9SPeter Brune Level: developer 780ab8d36c9SPeter Brune 781ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 782ab8d36c9SPeter Brune 783ab8d36c9SPeter Brune .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale() 784ab8d36c9SPeter Brune @*/ 785ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 786ab8d36c9SPeter Brune { 787ab8d36c9SPeter Brune SNES_FAS *fas; 7885fd66863SKarl Rupp 789ab8d36c9SPeter Brune PetscFunctionBegin; 790ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 791ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 792ab8d36c9SPeter Brune *vec = fas->rscale; 793ab8d36c9SPeter Brune PetscFunctionReturn(0); 794ab8d36c9SPeter Brune } 795ab8d36c9SPeter Brune 796ab8d36c9SPeter Brune #undef __FUNCT__ 797ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleIsFine" 798ab8d36c9SPeter Brune /*@ 799ab8d36c9SPeter Brune SNESFASCycleIsFine - Determines if a given cycle is the fine level. 800ab8d36c9SPeter Brune 801ab8d36c9SPeter Brune Logically Collective on SNES 802ab8d36c9SPeter Brune 803ab8d36c9SPeter Brune Input Parameters: 804ab8d36c9SPeter Brune . snes - the FAS context 805ab8d36c9SPeter Brune 806ab8d36c9SPeter Brune Output Parameters: 807ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not 808ab8d36c9SPeter Brune 809ab8d36c9SPeter Brune Level: advanced 810ab8d36c9SPeter Brune 811ab8d36c9SPeter Brune .keywords: SNES, FAS 812ab8d36c9SPeter Brune 813ab8d36c9SPeter Brune .seealso: SNESFASSetLevels() 814ab8d36c9SPeter Brune @*/ 815ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 816ab8d36c9SPeter Brune { 817ab8d36c9SPeter Brune SNES_FAS *fas; 8185fd66863SKarl Rupp 819ab8d36c9SPeter Brune PetscFunctionBegin; 820ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 821ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 8221aa26658SKarl Rupp if (fas->level == fas->levels - 1) *flg = PETSC_TRUE; 8231aa26658SKarl Rupp else *flg = PETSC_FALSE; 824ab8d36c9SPeter Brune PetscFunctionReturn(0); 825ab8d36c9SPeter Brune } 826ab8d36c9SPeter Brune 827ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */ 828ab8d36c9SPeter Brune 829ab8d36c9SPeter Brune #undef __FUNCT__ 830ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 831ab8d36c9SPeter Brune /*@ 832ab8d36c9SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 833ab8d36c9SPeter Brune interpolation from l-1 to the lth level 834ab8d36c9SPeter Brune 835ab8d36c9SPeter Brune Input Parameters: 836ab8d36c9SPeter Brune + snes - the multigrid context 837ab8d36c9SPeter Brune . mat - the interpolation operator 838ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 839ab8d36c9SPeter Brune 840ab8d36c9SPeter Brune Level: advanced 841ab8d36c9SPeter Brune 842ab8d36c9SPeter Brune Notes: 843ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction 844ab8d36c9SPeter Brune for the same level. 845ab8d36c9SPeter Brune 846ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 847ab8d36c9SPeter Brune out from the matrix size which one it is. 848ab8d36c9SPeter Brune 849ab8d36c9SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 850ab8d36c9SPeter Brune 851ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 852ab8d36c9SPeter Brune @*/ 8530adebc6cSBarry Smith PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) 8540adebc6cSBarry Smith { 85522d28d08SBarry Smith SNES_FAS *fas; 856ab8d36c9SPeter Brune PetscErrorCode ierr; 857ab8d36c9SPeter Brune SNES levelsnes; 85822d28d08SBarry Smith 859ab8d36c9SPeter Brune PetscFunctionBegin; 860ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 861ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 862ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 863ab8d36c9SPeter Brune ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 8641aa26658SKarl Rupp 865ab8d36c9SPeter Brune fas->interpolate = mat; 866ab8d36c9SPeter Brune PetscFunctionReturn(0); 867ab8d36c9SPeter Brune } 868ab8d36c9SPeter Brune 869ab8d36c9SPeter Brune #undef __FUNCT__ 870ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInterpolation" 871ab8d36c9SPeter Brune /*@ 872ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the 873ab8d36c9SPeter Brune interpolation from l-1 to the lth level 874ab8d36c9SPeter Brune 875ab8d36c9SPeter Brune Input Parameters: 876ab8d36c9SPeter Brune + snes - the multigrid context 877ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 878ab8d36c9SPeter Brune 879ab8d36c9SPeter Brune Output Parameters: 880ab8d36c9SPeter Brune . mat - the interpolation operator 881ab8d36c9SPeter Brune 882ab8d36c9SPeter Brune Level: advanced 883ab8d36c9SPeter Brune 884ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, interpolate, level 885ab8d36c9SPeter Brune 886ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale() 887ab8d36c9SPeter Brune @*/ 88822d28d08SBarry Smith PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) 88922d28d08SBarry Smith { 89022d28d08SBarry Smith SNES_FAS *fas; 891ab8d36c9SPeter Brune PetscErrorCode ierr; 892ab8d36c9SPeter Brune SNES levelsnes; 89322d28d08SBarry Smith 894ab8d36c9SPeter Brune PetscFunctionBegin; 895ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 896ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 897ab8d36c9SPeter Brune *mat = fas->interpolate; 898ab8d36c9SPeter Brune PetscFunctionReturn(0); 899ab8d36c9SPeter Brune } 900ab8d36c9SPeter Brune 901ab8d36c9SPeter Brune #undef __FUNCT__ 902ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 903ab8d36c9SPeter Brune /*@ 904ab8d36c9SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 905ab8d36c9SPeter Brune from level l to l-1. 906ab8d36c9SPeter Brune 907ab8d36c9SPeter Brune Input Parameters: 908ab8d36c9SPeter Brune + snes - the multigrid context 909ab8d36c9SPeter Brune . mat - the restriction matrix 910ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 911ab8d36c9SPeter Brune 912ab8d36c9SPeter Brune Level: advanced 913ab8d36c9SPeter Brune 914ab8d36c9SPeter Brune Notes: 915ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation 916ab8d36c9SPeter Brune for the same level. 917ab8d36c9SPeter Brune 918ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 919ab8d36c9SPeter Brune out from the matrix size which one it is. 920ab8d36c9SPeter Brune 921ab8d36c9SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 922ab8d36c9SPeter Brune is used. 923ab8d36c9SPeter Brune 924ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 925ab8d36c9SPeter Brune 926ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 927ab8d36c9SPeter Brune @*/ 92822d28d08SBarry Smith PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) 92922d28d08SBarry Smith { 93022d28d08SBarry Smith SNES_FAS *fas; 931ab8d36c9SPeter Brune PetscErrorCode ierr; 932ab8d36c9SPeter Brune SNES levelsnes; 93322d28d08SBarry Smith 934ab8d36c9SPeter Brune PetscFunctionBegin; 935ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 936ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 937ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 938ab8d36c9SPeter Brune ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 9391aa26658SKarl Rupp 940ab8d36c9SPeter Brune fas->restrct = mat; 941ab8d36c9SPeter Brune PetscFunctionReturn(0); 942ab8d36c9SPeter Brune } 943ab8d36c9SPeter Brune 944ab8d36c9SPeter Brune #undef __FUNCT__ 945ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetRestriction" 946ab8d36c9SPeter Brune /*@ 947ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the 948ab8d36c9SPeter Brune restriction from l to the l-1th level 949ab8d36c9SPeter Brune 950ab8d36c9SPeter Brune Input Parameters: 951ab8d36c9SPeter Brune + snes - the multigrid context 952ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 953ab8d36c9SPeter Brune 954ab8d36c9SPeter Brune Output Parameters: 955ab8d36c9SPeter Brune . mat - the interpolation operator 956ab8d36c9SPeter Brune 957ab8d36c9SPeter Brune Level: advanced 958ab8d36c9SPeter Brune 959ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, restrict, level 960ab8d36c9SPeter Brune 961ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale() 962ab8d36c9SPeter Brune @*/ 96322d28d08SBarry Smith PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) 96422d28d08SBarry Smith { 96522d28d08SBarry Smith SNES_FAS *fas; 966ab8d36c9SPeter Brune PetscErrorCode ierr; 967ab8d36c9SPeter Brune SNES levelsnes; 96822d28d08SBarry Smith 969ab8d36c9SPeter Brune PetscFunctionBegin; 970ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 971ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 972ab8d36c9SPeter Brune *mat = fas->restrct; 973ab8d36c9SPeter Brune PetscFunctionReturn(0); 974ab8d36c9SPeter Brune } 975ab8d36c9SPeter Brune 976ab8d36c9SPeter Brune 977ab8d36c9SPeter Brune #undef __FUNCT__ 978ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInjection" 979ab8d36c9SPeter Brune /*@ 980ab8d36c9SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 981ab8d36c9SPeter Brune from level l to l-1. 982ab8d36c9SPeter Brune 983ab8d36c9SPeter Brune Input Parameters: 984ab8d36c9SPeter Brune + snes - the multigrid context 985ab8d36c9SPeter Brune . mat - the restriction matrix 986ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 987ab8d36c9SPeter Brune 988ab8d36c9SPeter Brune Level: advanced 989ab8d36c9SPeter Brune 990ab8d36c9SPeter Brune Notes: 991ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to 992ab8d36c9SPeter Brune project the solution instead. 993ab8d36c9SPeter Brune 994ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 995ab8d36c9SPeter Brune 996ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 997ab8d36c9SPeter Brune @*/ 99822d28d08SBarry Smith PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) 99922d28d08SBarry Smith { 100022d28d08SBarry Smith SNES_FAS *fas; 1001ab8d36c9SPeter Brune PetscErrorCode ierr; 1002ab8d36c9SPeter Brune SNES levelsnes; 100322d28d08SBarry Smith 1004ab8d36c9SPeter Brune PetscFunctionBegin; 1005ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1006ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1007ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1008ab8d36c9SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 10091aa26658SKarl Rupp 1010ab8d36c9SPeter Brune fas->inject = mat; 1011ab8d36c9SPeter Brune PetscFunctionReturn(0); 1012ab8d36c9SPeter Brune } 1013ab8d36c9SPeter Brune 1014ab8d36c9SPeter Brune 1015ab8d36c9SPeter Brune #undef __FUNCT__ 1016ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInjection" 1017ab8d36c9SPeter Brune /*@ 1018ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the 1019ab8d36c9SPeter Brune injection from l-1 to the lth level 1020ab8d36c9SPeter Brune 1021ab8d36c9SPeter Brune Input Parameters: 1022ab8d36c9SPeter Brune + snes - the multigrid context 1023ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 1024ab8d36c9SPeter Brune 1025ab8d36c9SPeter Brune Output Parameters: 1026ab8d36c9SPeter Brune . mat - the injection operator 1027ab8d36c9SPeter Brune 1028ab8d36c9SPeter Brune Level: advanced 1029ab8d36c9SPeter Brune 1030ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, restrict, level 1031ab8d36c9SPeter Brune 1032ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale() 1033ab8d36c9SPeter Brune @*/ 103422d28d08SBarry Smith PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) 103522d28d08SBarry Smith { 103622d28d08SBarry Smith SNES_FAS *fas; 1037ab8d36c9SPeter Brune PetscErrorCode ierr; 1038ab8d36c9SPeter Brune SNES levelsnes; 103922d28d08SBarry Smith 1040ab8d36c9SPeter Brune PetscFunctionBegin; 1041ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1042ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1043ab8d36c9SPeter Brune *mat = fas->inject; 1044ab8d36c9SPeter Brune PetscFunctionReturn(0); 1045ab8d36c9SPeter Brune } 1046ab8d36c9SPeter Brune 1047ab8d36c9SPeter Brune #undef __FUNCT__ 1048ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 1049ab8d36c9SPeter Brune /*@ 1050ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 1051ab8d36c9SPeter Brune operator from level l to l-1. 1052ab8d36c9SPeter Brune 1053ab8d36c9SPeter Brune Input Parameters: 1054ab8d36c9SPeter Brune + snes - the multigrid context 1055ab8d36c9SPeter Brune . rscale - the restriction scaling 1056ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 1057ab8d36c9SPeter Brune 1058ab8d36c9SPeter Brune Level: advanced 1059ab8d36c9SPeter Brune 1060ab8d36c9SPeter Brune Notes: 1061ab8d36c9SPeter Brune This is only used in the case that the injection is not set. 1062ab8d36c9SPeter Brune 1063ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 1064ab8d36c9SPeter Brune 1065ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1066ab8d36c9SPeter Brune @*/ 106722d28d08SBarry Smith PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) 106822d28d08SBarry Smith { 1069ab8d36c9SPeter Brune SNES_FAS *fas; 1070ab8d36c9SPeter Brune PetscErrorCode ierr; 1071ab8d36c9SPeter Brune SNES levelsnes; 107222d28d08SBarry Smith 1073ab8d36c9SPeter Brune PetscFunctionBegin; 1074ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1075ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1076ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 1077ab8d36c9SPeter Brune ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 10781aa26658SKarl Rupp 1079ab8d36c9SPeter Brune fas->rscale = rscale; 1080ab8d36c9SPeter Brune PetscFunctionReturn(0); 1081ab8d36c9SPeter Brune } 1082ab8d36c9SPeter Brune 1083ab8d36c9SPeter Brune #undef __FUNCT__ 1084ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmoother" 1085ab8d36c9SPeter Brune /*@ 1086ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level. 1087ab8d36c9SPeter Brune 1088ab8d36c9SPeter Brune Input Parameters: 1089ab8d36c9SPeter Brune + snes - the multigrid context 1090ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1091ab8d36c9SPeter Brune 1092ab8d36c9SPeter Brune Output Parameters: 1093ab8d36c9SPeter Brune smooth - the smoother 1094ab8d36c9SPeter Brune 1095ab8d36c9SPeter Brune Level: advanced 1096ab8d36c9SPeter Brune 1097ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1098ab8d36c9SPeter Brune 1099ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1100ab8d36c9SPeter Brune @*/ 110122d28d08SBarry Smith PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) 110222d28d08SBarry Smith { 1103ab8d36c9SPeter Brune SNES_FAS *fas; 1104ab8d36c9SPeter Brune PetscErrorCode ierr; 1105ab8d36c9SPeter Brune SNES levelsnes; 110622d28d08SBarry Smith 1107ab8d36c9SPeter Brune PetscFunctionBegin; 1108ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1109ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1110ab8d36c9SPeter Brune if (!fas->smoothd) { 11113de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1112ab8d36c9SPeter Brune } 1113ab8d36c9SPeter Brune *smooth = fas->smoothd; 1114ab8d36c9SPeter Brune PetscFunctionReturn(0); 1115ab8d36c9SPeter Brune } 1116ab8d36c9SPeter Brune 1117ab8d36c9SPeter Brune #undef __FUNCT__ 1118ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherDown" 1119ab8d36c9SPeter Brune /*@ 1120ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level. 1121ab8d36c9SPeter Brune 1122ab8d36c9SPeter Brune Input Parameters: 1123ab8d36c9SPeter Brune + snes - the multigrid context 1124ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1125ab8d36c9SPeter Brune 1126ab8d36c9SPeter Brune Output Parameters: 1127ab8d36c9SPeter Brune smooth - the smoother 1128ab8d36c9SPeter Brune 1129ab8d36c9SPeter Brune Level: advanced 1130ab8d36c9SPeter Brune 1131ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1132ab8d36c9SPeter Brune 1133ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1134ab8d36c9SPeter Brune @*/ 113522d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) 113622d28d08SBarry Smith { 1137ab8d36c9SPeter Brune SNES_FAS *fas; 1138ab8d36c9SPeter Brune PetscErrorCode ierr; 1139ab8d36c9SPeter Brune SNES levelsnes; 114022d28d08SBarry Smith 1141ab8d36c9SPeter Brune PetscFunctionBegin; 1142ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1143ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1144ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1145ab8d36c9SPeter Brune if (!fas->smoothd) { 11463de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1147ab8d36c9SPeter Brune } 1148ab8d36c9SPeter Brune if (!fas->smoothu) { 11493de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr); 1150ab8d36c9SPeter Brune } 1151ab8d36c9SPeter Brune *smooth = fas->smoothd; 1152ab8d36c9SPeter Brune PetscFunctionReturn(0); 1153ab8d36c9SPeter Brune } 1154ab8d36c9SPeter Brune 1155ab8d36c9SPeter Brune #undef __FUNCT__ 1156ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherUp" 1157ab8d36c9SPeter Brune /*@ 1158ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level. 1159ab8d36c9SPeter Brune 1160ab8d36c9SPeter Brune Input Parameters: 1161ab8d36c9SPeter Brune + snes - the multigrid context 1162ab8d36c9SPeter Brune - level - the level (0 is coarsest) 1163ab8d36c9SPeter Brune 1164ab8d36c9SPeter Brune Output Parameters: 1165ab8d36c9SPeter Brune smooth - the smoother 1166ab8d36c9SPeter Brune 1167ab8d36c9SPeter Brune Level: advanced 1168ab8d36c9SPeter Brune 1169ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1170ab8d36c9SPeter Brune 1171ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1172ab8d36c9SPeter Brune @*/ 117322d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) 117422d28d08SBarry Smith { 1175ab8d36c9SPeter Brune SNES_FAS *fas; 1176ab8d36c9SPeter Brune PetscErrorCode ierr; 1177ab8d36c9SPeter Brune SNES levelsnes; 117822d28d08SBarry Smith 1179ab8d36c9SPeter Brune PetscFunctionBegin; 1180ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1181ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1182ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1183ab8d36c9SPeter Brune if (!fas->smoothd) { 11843de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1185ab8d36c9SPeter Brune } 1186ab8d36c9SPeter Brune if (!fas->smoothu) { 11873de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr); 1188ab8d36c9SPeter Brune } 1189ab8d36c9SPeter Brune *smooth = fas->smoothu; 1190ab8d36c9SPeter Brune PetscFunctionReturn(0); 1191ab8d36c9SPeter Brune } 1192d6ad1212SPeter Brune 1193d6ad1212SPeter Brune #undef __FUNCT__ 1194d6ad1212SPeter Brune #define __FUNCT__ "SNESFASGetCoarseSolve" 1195d6ad1212SPeter Brune /*@ 1196d6ad1212SPeter Brune SNESFASGetCoarseSolve - Gets the coarsest solver. 1197d6ad1212SPeter Brune 1198d6ad1212SPeter Brune Input Parameters: 1199a3a80b83SMatthew G. Knepley . snes - the multigrid context 1200d6ad1212SPeter Brune 1201d6ad1212SPeter Brune Output Parameters: 1202a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver 1203d6ad1212SPeter Brune 1204d6ad1212SPeter Brune Level: advanced 1205d6ad1212SPeter Brune 1206d6ad1212SPeter Brune .keywords: FAS, MG, get, multigrid, solver, coarse 1207d6ad1212SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1208d6ad1212SPeter Brune @*/ 1209a3a80b83SMatthew G. Knepley PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse) 121022d28d08SBarry Smith { 1211d6ad1212SPeter Brune SNES_FAS *fas; 1212d6ad1212SPeter Brune PetscErrorCode ierr; 1213d6ad1212SPeter Brune SNES levelsnes; 121422d28d08SBarry Smith 1215d6ad1212SPeter Brune PetscFunctionBegin; 1216d6ad1212SPeter Brune ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr); 1217d6ad1212SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1218d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1219d6ad1212SPeter Brune if (!fas->smoothd) { 12203de4d1dbSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr); 1221d6ad1212SPeter Brune } 1222a3a80b83SMatthew G. Knepley *coarse = fas->smoothd; 1223d6ad1212SPeter Brune PetscFunctionReturn(0); 1224d6ad1212SPeter Brune } 1225928e959bSPeter Brune 1226928e959bSPeter Brune #undef __FUNCT__ 1227928e959bSPeter Brune #define __FUNCT__ "SNESFASFullSetDownSweep" 1228928e959bSPeter Brune /*@ 1229928e959bSPeter Brune SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS 1230928e959bSPeter Brune 1231928e959bSPeter Brune Logically Collective on SNES 1232928e959bSPeter Brune 1233928e959bSPeter Brune Input Parameters: 1234928e959bSPeter Brune + snes - the multigrid context 1235928e959bSPeter Brune - swp - whether to downsweep or not 1236928e959bSPeter Brune 1237928e959bSPeter Brune Options Database Key: 1238928e959bSPeter Brune . -snes_fas_full_downsweep - Sets number of pre-smoothing steps 1239928e959bSPeter Brune 1240928e959bSPeter Brune Level: advanced 1241928e959bSPeter Brune 1242928e959bSPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 1243928e959bSPeter Brune 1244928e959bSPeter Brune .seealso: SNESFASSetNumberSmoothUp() 1245928e959bSPeter Brune @*/ 1246928e959bSPeter Brune PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp) 1247928e959bSPeter Brune { 1248928e959bSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 1249928e959bSPeter Brune PetscErrorCode ierr = 0; 1250928e959bSPeter Brune 1251928e959bSPeter Brune PetscFunctionBegin; 1252928e959bSPeter Brune fas->full_downsweep = swp; 1253928e959bSPeter Brune if (fas->next) { 1254928e959bSPeter Brune ierr = SNESFASFullSetDownSweep(fas->next,swp);CHKERRQ(ierr); 1255928e959bSPeter Brune } 1256928e959bSPeter Brune PetscFunctionReturn(0); 1257928e959bSPeter Brune } 1258