1ab8d36c9SPeter Brune #include "../src/snes/impls/fas/fasimpls.h" /*I "petscsnesfas.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); 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); 1291aa26658SKarl Rupp 130ab8d36c9SPeter Brune ((SNES_FAS*)fas->next->data)->previous = prevsnes; 1311aa26658SKarl Rupp 132ab8d36c9SPeter Brune prevsnes = fas->next; 133ab8d36c9SPeter Brune fas = (SNES_FAS*)prevsnes->data; 134ab8d36c9SPeter Brune } 135ab8d36c9SPeter Brune } 136ab8d36c9SPeter Brune PetscFunctionReturn(0); 137ab8d36c9SPeter Brune } 138ab8d36c9SPeter Brune 139ab8d36c9SPeter Brune 140ab8d36c9SPeter Brune #undef __FUNCT__ 141ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 142ab8d36c9SPeter Brune /*@ 143ab8d36c9SPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 144ab8d36c9SPeter Brune 145ab8d36c9SPeter Brune Input Parameter: 146ab8d36c9SPeter Brune . snes - the nonlinear solver context 147ab8d36c9SPeter Brune 148ab8d36c9SPeter Brune Output parameter: 149ab8d36c9SPeter Brune . levels - the number of levels 150ab8d36c9SPeter Brune 151ab8d36c9SPeter Brune Level: advanced 152ab8d36c9SPeter Brune 153ab8d36c9SPeter Brune .keywords: MG, get, levels, multigrid 154ab8d36c9SPeter Brune 155ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 156ab8d36c9SPeter Brune @*/ 15722d28d08SBarry Smith PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels) 15822d28d08SBarry Smith { 159ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 1605fd66863SKarl Rupp 161ab8d36c9SPeter Brune PetscFunctionBegin; 162ab8d36c9SPeter Brune *levels = fas->levels; 163ab8d36c9SPeter Brune PetscFunctionReturn(0); 164ab8d36c9SPeter Brune } 165ab8d36c9SPeter Brune 166ab8d36c9SPeter Brune 167ab8d36c9SPeter Brune #undef __FUNCT__ 168ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetCycleSNES" 169ab8d36c9SPeter Brune /*@ 170ab8d36c9SPeter Brune SNESFASGetCycleSNES - Gets the SNES corresponding to a particular 171ab8d36c9SPeter Brune level of the FAS hierarchy. 172ab8d36c9SPeter Brune 173ab8d36c9SPeter Brune Input Parameters: 174ab8d36c9SPeter Brune + snes - the multigrid context 175ab8d36c9SPeter Brune level - the level to get 176ab8d36c9SPeter Brune - lsnes - whether to use the nonlinear smoother or not 177ab8d36c9SPeter Brune 178ab8d36c9SPeter Brune Level: advanced 179ab8d36c9SPeter Brune 180ab8d36c9SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 181ab8d36c9SPeter Brune 182ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 183ab8d36c9SPeter Brune @*/ 18422d28d08SBarry Smith PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes) 18522d28d08SBarry Smith { 186ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 187ab8d36c9SPeter Brune PetscInt i; 188ab8d36c9SPeter Brune 189ab8d36c9SPeter Brune PetscFunctionBegin; 190ce94432eSBarry 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); 191ce94432eSBarry 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); 192ab8d36c9SPeter Brune 193ab8d36c9SPeter Brune *lsnes = snes; 194ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) { 195ab8d36c9SPeter Brune *lsnes = fas->next; 196ab8d36c9SPeter Brune fas = (SNES_FAS*)(*lsnes)->data; 197ab8d36c9SPeter Brune } 198ce94432eSBarry Smith if (fas->level != level) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 199ab8d36c9SPeter Brune PetscFunctionReturn(0); 200ab8d36c9SPeter Brune } 201ab8d36c9SPeter Brune 202ab8d36c9SPeter Brune #undef __FUNCT__ 203ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 204ab8d36c9SPeter Brune /*@ 205ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 206ab8d36c9SPeter Brune use on all levels. 207ab8d36c9SPeter Brune 208ab8d36c9SPeter Brune Logically Collective on SNES 209ab8d36c9SPeter Brune 210ab8d36c9SPeter Brune Input Parameters: 211ab8d36c9SPeter Brune + snes - the multigrid context 212ab8d36c9SPeter Brune - n - the number of smoothing steps 213ab8d36c9SPeter Brune 214ab8d36c9SPeter Brune Options Database Key: 215ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 216ab8d36c9SPeter Brune 217ab8d36c9SPeter Brune Level: advanced 218ab8d36c9SPeter Brune 219ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 220ab8d36c9SPeter Brune 221ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 222ab8d36c9SPeter Brune @*/ 22322d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) 22422d28d08SBarry Smith { 225ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 22622d28d08SBarry Smith PetscErrorCode ierr; 22722d28d08SBarry Smith 228ab8d36c9SPeter Brune PetscFunctionBegin; 229ab8d36c9SPeter Brune fas->max_up_it = n; 230656ede7eSPeter Brune if (!fas->smoothu && fas->level != 0) { 23122d28d08SBarry Smith ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 232ab8d36c9SPeter Brune } 23322d28d08SBarry Smith if (fas->smoothu) { 23422d28d08SBarry Smith ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr); 23522d28d08SBarry Smith } 236ab8d36c9SPeter Brune if (fas->next) { 237ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 238ab8d36c9SPeter Brune } 239ab8d36c9SPeter Brune PetscFunctionReturn(0); 240ab8d36c9SPeter Brune } 241ab8d36c9SPeter Brune 242ab8d36c9SPeter Brune #undef __FUNCT__ 243ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 244ab8d36c9SPeter Brune /*@ 245ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 246ab8d36c9SPeter Brune use on all levels. 247ab8d36c9SPeter Brune 248ab8d36c9SPeter Brune Logically Collective on SNES 249ab8d36c9SPeter Brune 250ab8d36c9SPeter Brune Input Parameters: 251ab8d36c9SPeter Brune + snes - the multigrid context 252ab8d36c9SPeter Brune - n - the number of smoothing steps 253ab8d36c9SPeter Brune 254ab8d36c9SPeter Brune Options Database Key: 255ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 256ab8d36c9SPeter Brune 257ab8d36c9SPeter Brune Level: advanced 258ab8d36c9SPeter Brune 259ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 260ab8d36c9SPeter Brune 261ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 262ab8d36c9SPeter Brune @*/ 26322d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) 26422d28d08SBarry Smith { 265ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 266ab8d36c9SPeter Brune PetscErrorCode ierr = 0; 26722d28d08SBarry Smith 268ab8d36c9SPeter Brune PetscFunctionBegin; 269ab8d36c9SPeter Brune if (!fas->smoothd) { 27022d28d08SBarry Smith ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 271ab8d36c9SPeter Brune } 272ab8d36c9SPeter Brune ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr); 2731aa26658SKarl Rupp 274ab8d36c9SPeter Brune fas->max_down_it = n; 275ab8d36c9SPeter Brune if (fas->next) { 276ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 277ab8d36c9SPeter Brune } 278ab8d36c9SPeter Brune PetscFunctionReturn(0); 279ab8d36c9SPeter Brune } 280ab8d36c9SPeter Brune 281ab8d36c9SPeter Brune 282ab8d36c9SPeter Brune #undef __FUNCT__ 283*87f44e3fSPeter Brune #define __FUNCT__ "SNESFASSetContinuation" 284*87f44e3fSPeter Brune /*@ 285*87f44e3fSPeter Brune SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep 286*87f44e3fSPeter Brune 287*87f44e3fSPeter Brune Logically Collective on SNES 288*87f44e3fSPeter Brune 289*87f44e3fSPeter Brune Input Parameters: 290*87f44e3fSPeter Brune + snes - the multigrid context 291*87f44e3fSPeter Brune - n - the number of smoothing steps 292*87f44e3fSPeter Brune 293*87f44e3fSPeter Brune Options Database Key: 294*87f44e3fSPeter Brune . -snes_fas_continuation - sets continuation to true 295*87f44e3fSPeter Brune 296*87f44e3fSPeter Brune Level: advanced 297*87f44e3fSPeter Brune 298*87f44e3fSPeter Brune Notes: This sets the prefix on the upsweep smoothers to -fas_continuation 299*87f44e3fSPeter Brune 300*87f44e3fSPeter Brune .keywords: FAS, MG, smoother, continuation 301*87f44e3fSPeter Brune 302*87f44e3fSPeter Brune .seealso: SNESFAS 303*87f44e3fSPeter Brune @*/ 304*87f44e3fSPeter Brune PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation) 305*87f44e3fSPeter Brune { 306*87f44e3fSPeter Brune const char *optionsprefix; 307*87f44e3fSPeter Brune char tprefix[128]; 308*87f44e3fSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 309*87f44e3fSPeter Brune PetscErrorCode ierr = 0; 310*87f44e3fSPeter Brune 311*87f44e3fSPeter Brune PetscFunctionBegin; 312*87f44e3fSPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 313*87f44e3fSPeter Brune if (!fas->smoothu) { 314*87f44e3fSPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 315*87f44e3fSPeter Brune } 316*87f44e3fSPeter Brune sprintf(tprefix,"fas_levels_continuation_"); 317*87f44e3fSPeter Brune ierr = SNESSetOptionsPrefix(fas->smoothu, optionsprefix);CHKERRQ(ierr); 318*87f44e3fSPeter Brune ierr = SNESAppendOptionsPrefix(fas->smoothu, tprefix);CHKERRQ(ierr); 319*87f44e3fSPeter Brune ierr = SNESSetType(fas->smoothu,SNESNEWTONLS);CHKERRQ(ierr); 320*87f44e3fSPeter Brune ierr = SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100);CHKERRQ(ierr); 321*87f44e3fSPeter Brune fas->continuation = continuation; 322*87f44e3fSPeter Brune if (fas->next) { 323*87f44e3fSPeter Brune ierr = SNESFASSetContinuation(fas->next,continuation);CHKERRQ(ierr); 324*87f44e3fSPeter Brune } 325*87f44e3fSPeter Brune PetscFunctionReturn(0); 326*87f44e3fSPeter Brune } 327*87f44e3fSPeter Brune 328*87f44e3fSPeter Brune 329*87f44e3fSPeter Brune #undef __FUNCT__ 330ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetCycles" 331ab8d36c9SPeter Brune /*@ 332ab8d36c9SPeter Brune SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 333ab8d36c9SPeter Brune complicated cycling. 334ab8d36c9SPeter Brune 335ab8d36c9SPeter Brune Logically Collective on SNES 336ab8d36c9SPeter Brune 337ab8d36c9SPeter Brune Input Parameters: 338ab8d36c9SPeter Brune + snes - the multigrid context 339ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 340ab8d36c9SPeter Brune 341ab8d36c9SPeter Brune Options Database Key: 342ab8d36c9SPeter Brune $ -snes_fas_cycles 1 or 2 343ab8d36c9SPeter Brune 344ab8d36c9SPeter Brune Level: advanced 345ab8d36c9SPeter Brune 346ab8d36c9SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 347ab8d36c9SPeter Brune 348ab8d36c9SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 349ab8d36c9SPeter Brune @*/ 35022d28d08SBarry Smith PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) 35122d28d08SBarry Smith { 352ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 353ab8d36c9SPeter Brune PetscErrorCode ierr; 354ab8d36c9SPeter Brune PetscBool isFine; 35522d28d08SBarry Smith 356ab8d36c9SPeter Brune PetscFunctionBegin; 35722d28d08SBarry Smith ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 3581aa26658SKarl Rupp 359ab8d36c9SPeter Brune fas->n_cycles = cycles; 3601aa26658SKarl Rupp if (!isFine) { 361ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 3621aa26658SKarl Rupp } 363ab8d36c9SPeter Brune if (fas->next) { 364ab8d36c9SPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 365ab8d36c9SPeter Brune } 366ab8d36c9SPeter Brune PetscFunctionReturn(0); 367ab8d36c9SPeter Brune } 368ab8d36c9SPeter Brune 369c8c899caSPeter Brune 370c8c899caSPeter Brune #undef __FUNCT__ 371c8c899caSPeter Brune #define __FUNCT__ "SNESFASSetMonitor" 372c8c899caSPeter Brune /*@ 373c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring 374c8c899caSPeter Brune 375c8c899caSPeter Brune Logically Collective on SNES 376c8c899caSPeter Brune 377c8c899caSPeter Brune Input Parameters: 378c8c899caSPeter Brune + snes - the FAS context 379c8c899caSPeter Brune - flg - monitor or not 380c8c899caSPeter Brune 381c8c899caSPeter Brune Level: advanced 382c8c899caSPeter Brune 383c8c899caSPeter Brune .keywords: FAS, monitor 384c8c899caSPeter Brune 385c8c899caSPeter Brune .seealso: SNESFASSetCyclesOnLevel() 386c8c899caSPeter Brune @*/ 38722d28d08SBarry Smith PetscErrorCode SNESFASSetMonitor(SNES snes, PetscBool flg) 38822d28d08SBarry Smith { 389c8c899caSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 390c8c899caSPeter Brune PetscErrorCode ierr; 391c8c899caSPeter Brune PetscBool isFine; 392c8c899caSPeter Brune PetscInt i, levels = fas->levels; 393c8c899caSPeter Brune SNES levelsnes; 39422d28d08SBarry Smith 395c8c899caSPeter Brune PetscFunctionBegin; 396c8c899caSPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 397c8c899caSPeter Brune if (isFine) { 398c8c899caSPeter Brune for (i = 0; i < levels; i++) { 39922d28d08SBarry Smith ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 400c8c899caSPeter Brune fas = (SNES_FAS*)levelsnes->data; 401c8c899caSPeter Brune if (flg) { 402ce94432eSBarry Smith fas->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)levelsnes));CHKERRQ(ierr); 403c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */ 404c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 4050298fd71SBarry Smith ierr = SNESMonitorSet(levelsnes,SNESMonitorDefault,NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 4061aa26658SKarl Rupp } else if (i != fas->levels - 1) { 407c8c899caSPeter Brune /* unset the monitors on the coarse levels */ 408c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 409c8c899caSPeter Brune } 410c8c899caSPeter Brune } 411c8c899caSPeter Brune } 412c8c899caSPeter Brune PetscFunctionReturn(0); 413c8c899caSPeter Brune } 414c8c899caSPeter Brune 415ab8d36c9SPeter Brune #undef __FUNCT__ 4160dd27c6cSPeter Brune #define __FUNCT__ "SNESFASSetLog" 4170dd27c6cSPeter Brune /*@ 4180dd27c6cSPeter Brune SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels 4190dd27c6cSPeter Brune 4200dd27c6cSPeter Brune Logically Collective on SNES 4210dd27c6cSPeter Brune 4220dd27c6cSPeter Brune Input Parameters: 4230dd27c6cSPeter Brune + snes - the FAS context 4240dd27c6cSPeter Brune - flg - monitor or not 4250dd27c6cSPeter Brune 4260dd27c6cSPeter Brune Level: advanced 4270dd27c6cSPeter Brune 4280dd27c6cSPeter Brune .keywords: FAS, logging 4290dd27c6cSPeter Brune 4300dd27c6cSPeter Brune .seealso: SNESFASSetMonitor() 4310dd27c6cSPeter Brune @*/ 4320dd27c6cSPeter Brune PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) 4330dd27c6cSPeter Brune { 4340dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 4350dd27c6cSPeter Brune PetscErrorCode ierr; 4360dd27c6cSPeter Brune PetscBool isFine; 4370dd27c6cSPeter Brune PetscInt i, levels = fas->levels; 4380dd27c6cSPeter Brune SNES levelsnes; 4390dd27c6cSPeter Brune char eventname[128]; 4400dd27c6cSPeter Brune 4410dd27c6cSPeter Brune PetscFunctionBegin; 4420dd27c6cSPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 4430dd27c6cSPeter Brune if (isFine) { 4440dd27c6cSPeter Brune for (i = 0; i < levels; i++) { 4450dd27c6cSPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 4460dd27c6cSPeter Brune fas = (SNES_FAS*)levelsnes->data; 4470dd27c6cSPeter Brune if (flg) { 4480dd27c6cSPeter Brune sprintf(eventname,"FASSetup %d",(int)i); 4490dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr); 4500dd27c6cSPeter Brune sprintf(eventname,"FASSmooth %d",(int)i); 4510dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr); 4520dd27c6cSPeter Brune sprintf(eventname,"FASResid %d",(int)i); 4530dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr); 4540dd27c6cSPeter Brune if (i) { 4550dd27c6cSPeter Brune sprintf(eventname,"FASInterp %d",(int)i); 4560dd27c6cSPeter Brune ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr); 4570dd27c6cSPeter Brune } 4580dd27c6cSPeter Brune } else { 4590298fd71SBarry Smith fas->eventsmoothsetup = 0; 4600298fd71SBarry Smith fas->eventsmoothsolve = 0; 4610298fd71SBarry Smith fas->eventresidual = 0; 4620298fd71SBarry Smith fas->eventinterprestrict = 0; 4630dd27c6cSPeter Brune } 4640dd27c6cSPeter Brune } 4650dd27c6cSPeter Brune } 4660dd27c6cSPeter Brune PetscFunctionReturn(0); 4670dd27c6cSPeter Brune } 4680dd27c6cSPeter Brune 4690dd27c6cSPeter Brune #undef __FUNCT__ 470ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleCreateSmoother_Private" 471ab8d36c9SPeter Brune /* 472ab8d36c9SPeter Brune Creates the default smoother type. 473ab8d36c9SPeter Brune 47404d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level. 475ab8d36c9SPeter Brune 476ab8d36c9SPeter Brune */ 47722d28d08SBarry Smith PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) 47822d28d08SBarry Smith { 479ab8d36c9SPeter Brune SNES_FAS *fas; 480ab8d36c9SPeter Brune const char *optionsprefix; 481ab8d36c9SPeter Brune char tprefix[128]; 482ab8d36c9SPeter Brune PetscErrorCode ierr; 483ab8d36c9SPeter Brune SNES nsmooth; 48422d28d08SBarry Smith 485ab8d36c9SPeter Brune PetscFunctionBegin; 486ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 487ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 488ab8d36c9SPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 489ab8d36c9SPeter Brune /* create the default smoother */ 490ce94432eSBarry Smith ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);CHKERRQ(ierr); 491ab8d36c9SPeter Brune if (fas->level == 0) { 492ab8d36c9SPeter Brune sprintf(tprefix,"fas_coarse_"); 493ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 494ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 49504d7464bSBarry Smith ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr); 496e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr); 497ab8d36c9SPeter Brune } else { 498ab8d36c9SPeter Brune sprintf(tprefix,"fas_levels_%d_",(int)fas->level); 499ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 500ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 501ab8d36c9SPeter Brune ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr); 502e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr); 503ab8d36c9SPeter Brune } 504ab8d36c9SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 5053bb1ff40SBarry Smith ierr = PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);CHKERRQ(ierr); 506f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr); 507ab8d36c9SPeter Brune *smooth = nsmooth; 508ab8d36c9SPeter Brune PetscFunctionReturn(0); 509ab8d36c9SPeter Brune } 510ab8d36c9SPeter Brune 511ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */ 512ab8d36c9SPeter Brune 513ab8d36c9SPeter Brune #undef __FUNCT__ 514ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleSetCycles" 515ab8d36c9SPeter Brune /*@ 516ab8d36c9SPeter Brune SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 517ab8d36c9SPeter Brune 518ab8d36c9SPeter Brune Logically Collective on SNES 519ab8d36c9SPeter Brune 520ab8d36c9SPeter Brune Input Parameters: 521ab8d36c9SPeter Brune + snes - the multigrid context 522ab8d36c9SPeter Brune . level - the level to set the number of cycles on 523ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 524ab8d36c9SPeter Brune 525ab8d36c9SPeter Brune Level: advanced 526ab8d36c9SPeter Brune 527ab8d36c9SPeter Brune .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid 528ab8d36c9SPeter Brune 529ab8d36c9SPeter Brune .seealso: SNESFASSetCycles() 530ab8d36c9SPeter Brune @*/ 53122d28d08SBarry Smith PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) 53222d28d08SBarry Smith { 533ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 534ab8d36c9SPeter Brune PetscErrorCode ierr; 53522d28d08SBarry Smith 536ab8d36c9SPeter Brune PetscFunctionBegin; 537ab8d36c9SPeter Brune fas->n_cycles = cycles; 538ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 539ab8d36c9SPeter Brune PetscFunctionReturn(0); 540ab8d36c9SPeter Brune } 541ab8d36c9SPeter Brune 542ab8d36c9SPeter Brune 543ab8d36c9SPeter Brune #undef __FUNCT__ 544ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmoother" 545ab8d36c9SPeter Brune /*@ 546ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 547ab8d36c9SPeter Brune 548ab8d36c9SPeter Brune Logically Collective on SNES 549ab8d36c9SPeter Brune 550ab8d36c9SPeter Brune Input Parameters: 551ab8d36c9SPeter Brune . snes - the multigrid context 552ab8d36c9SPeter Brune 553ab8d36c9SPeter Brune Output Parameters: 554ab8d36c9SPeter Brune . smooth - the smoother 555ab8d36c9SPeter Brune 556ab8d36c9SPeter Brune Level: advanced 557ab8d36c9SPeter Brune 558ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 559ab8d36c9SPeter Brune 560ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown() 561ab8d36c9SPeter Brune @*/ 562ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 563ab8d36c9SPeter Brune { 564ab8d36c9SPeter Brune SNES_FAS *fas; 5655fd66863SKarl Rupp 566ab8d36c9SPeter Brune PetscFunctionBegin; 567ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 568ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 569ab8d36c9SPeter Brune *smooth = fas->smoothd; 570ab8d36c9SPeter Brune PetscFunctionReturn(0); 571ab8d36c9SPeter Brune } 572ab8d36c9SPeter Brune #undef __FUNCT__ 573ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherUp" 574ab8d36c9SPeter Brune /*@ 575ab8d36c9SPeter Brune SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 576ab8d36c9SPeter Brune 577ab8d36c9SPeter Brune Logically Collective on SNES 578ab8d36c9SPeter Brune 579ab8d36c9SPeter Brune Input Parameters: 580ab8d36c9SPeter Brune . snes - the multigrid context 581ab8d36c9SPeter Brune 582ab8d36c9SPeter Brune Output Parameters: 583ab8d36c9SPeter Brune . smoothu - the smoother 584ab8d36c9SPeter Brune 585ab8d36c9SPeter Brune Notes: 586ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent 587ab8d36c9SPeter Brune default behavior in the process of the solve. 588ab8d36c9SPeter Brune 589ab8d36c9SPeter Brune Level: advanced 590ab8d36c9SPeter Brune 591ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 592ab8d36c9SPeter Brune 593ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown() 594ab8d36c9SPeter Brune @*/ 595ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 596ab8d36c9SPeter Brune { 597ab8d36c9SPeter Brune SNES_FAS *fas; 5985fd66863SKarl Rupp 599ab8d36c9SPeter Brune PetscFunctionBegin; 600ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 601ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 6021aa26658SKarl Rupp if (!fas->smoothu) *smoothu = fas->smoothd; 6031aa26658SKarl Rupp else *smoothu = fas->smoothu; 604ab8d36c9SPeter Brune PetscFunctionReturn(0); 605ab8d36c9SPeter Brune } 606ab8d36c9SPeter Brune 607ab8d36c9SPeter Brune #undef __FUNCT__ 608ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherDown" 609ab8d36c9SPeter Brune /*@ 610ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 611ab8d36c9SPeter Brune 612ab8d36c9SPeter Brune Logically Collective on SNES 613ab8d36c9SPeter Brune 614ab8d36c9SPeter Brune Input Parameters: 615ab8d36c9SPeter Brune . snes - the multigrid context 616ab8d36c9SPeter Brune 617ab8d36c9SPeter Brune Output Parameters: 618ab8d36c9SPeter Brune . smoothd - the smoother 619ab8d36c9SPeter Brune 620ab8d36c9SPeter Brune Level: advanced 621ab8d36c9SPeter Brune 622ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 623ab8d36c9SPeter Brune 624ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 625ab8d36c9SPeter Brune @*/ 626ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 627ab8d36c9SPeter Brune { 628ab8d36c9SPeter Brune SNES_FAS *fas; 6295fd66863SKarl Rupp 630ab8d36c9SPeter Brune PetscFunctionBegin; 631ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 632ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 633ab8d36c9SPeter Brune *smoothd = fas->smoothd; 634ab8d36c9SPeter Brune PetscFunctionReturn(0); 635ab8d36c9SPeter Brune } 636ab8d36c9SPeter Brune 637ab8d36c9SPeter Brune 638ab8d36c9SPeter Brune #undef __FUNCT__ 639ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetCorrection" 640ab8d36c9SPeter Brune /*@ 641ab8d36c9SPeter Brune SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 642ab8d36c9SPeter Brune 643ab8d36c9SPeter Brune Logically Collective on SNES 644ab8d36c9SPeter Brune 645ab8d36c9SPeter Brune Input Parameters: 646ab8d36c9SPeter Brune . snes - the multigrid context 647ab8d36c9SPeter Brune 648ab8d36c9SPeter Brune Output Parameters: 649ab8d36c9SPeter Brune . correction - the coarse correction on this level 650ab8d36c9SPeter Brune 651ab8d36c9SPeter Brune Notes: 6520298fd71SBarry Smith Returns NULL on the coarsest level. 653ab8d36c9SPeter Brune 654ab8d36c9SPeter Brune Level: advanced 655ab8d36c9SPeter Brune 656ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 657ab8d36c9SPeter Brune 658ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 659ab8d36c9SPeter Brune @*/ 660ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 661ab8d36c9SPeter Brune { 662ab8d36c9SPeter Brune SNES_FAS *fas; 6635fd66863SKarl Rupp 664ab8d36c9SPeter Brune PetscFunctionBegin; 665ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 666ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 667ab8d36c9SPeter Brune *correction = fas->next; 668ab8d36c9SPeter Brune PetscFunctionReturn(0); 669ab8d36c9SPeter Brune } 670ab8d36c9SPeter Brune 671ab8d36c9SPeter Brune #undef __FUNCT__ 672ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInterpolation" 673ab8d36c9SPeter Brune /*@ 674ab8d36c9SPeter Brune SNESFASCycleGetInterpolation - Gets the interpolation on this level 675ab8d36c9SPeter Brune 676ab8d36c9SPeter Brune Logically Collective on SNES 677ab8d36c9SPeter Brune 678ab8d36c9SPeter Brune Input Parameters: 679ab8d36c9SPeter Brune . snes - the multigrid context 680ab8d36c9SPeter Brune 681ab8d36c9SPeter Brune Output Parameters: 682ab8d36c9SPeter Brune . mat - the interpolation operator on this level 683ab8d36c9SPeter Brune 684ab8d36c9SPeter Brune Level: developer 685ab8d36c9SPeter Brune 686ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 687ab8d36c9SPeter Brune 688ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 689ab8d36c9SPeter Brune @*/ 690ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 691ab8d36c9SPeter Brune { 692ab8d36c9SPeter Brune SNES_FAS *fas; 6935fd66863SKarl Rupp 694ab8d36c9SPeter Brune PetscFunctionBegin; 695ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 696ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 697ab8d36c9SPeter Brune *mat = fas->interpolate; 698ab8d36c9SPeter Brune PetscFunctionReturn(0); 699ab8d36c9SPeter Brune } 700ab8d36c9SPeter Brune 701ab8d36c9SPeter Brune 702ab8d36c9SPeter Brune #undef __FUNCT__ 703ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRestriction" 704ab8d36c9SPeter Brune /*@ 705ab8d36c9SPeter Brune SNESFASCycleGetRestriction - Gets the restriction on this level 706ab8d36c9SPeter Brune 707ab8d36c9SPeter Brune Logically Collective on SNES 708ab8d36c9SPeter Brune 709ab8d36c9SPeter Brune Input Parameters: 710ab8d36c9SPeter Brune . snes - the multigrid context 711ab8d36c9SPeter Brune 712ab8d36c9SPeter Brune Output Parameters: 713ab8d36c9SPeter Brune . mat - the restriction operator on this level 714ab8d36c9SPeter Brune 715ab8d36c9SPeter Brune Level: developer 716ab8d36c9SPeter Brune 717ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 718ab8d36c9SPeter Brune 719ab8d36c9SPeter Brune .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation() 720ab8d36c9SPeter Brune @*/ 721ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 722ab8d36c9SPeter Brune { 723ab8d36c9SPeter Brune SNES_FAS *fas; 7245fd66863SKarl Rupp 725ab8d36c9SPeter Brune PetscFunctionBegin; 726ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 727ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 728ab8d36c9SPeter Brune *mat = fas->restrct; 729ab8d36c9SPeter Brune PetscFunctionReturn(0); 730ab8d36c9SPeter Brune } 731ab8d36c9SPeter Brune 732ab8d36c9SPeter Brune 733ab8d36c9SPeter Brune #undef __FUNCT__ 734ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInjection" 735ab8d36c9SPeter Brune /*@ 736ab8d36c9SPeter Brune SNESFASCycleGetInjection - Gets the injection on this level 737ab8d36c9SPeter Brune 738ab8d36c9SPeter Brune Logically Collective on SNES 739ab8d36c9SPeter Brune 740ab8d36c9SPeter Brune Input Parameters: 741ab8d36c9SPeter Brune . snes - the multigrid context 742ab8d36c9SPeter Brune 743ab8d36c9SPeter Brune Output Parameters: 744ab8d36c9SPeter Brune . mat - the restriction operator on this level 745ab8d36c9SPeter Brune 746ab8d36c9SPeter Brune Level: developer 747ab8d36c9SPeter Brune 748ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 749ab8d36c9SPeter Brune 750ab8d36c9SPeter Brune .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction() 751ab8d36c9SPeter Brune @*/ 752ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 753ab8d36c9SPeter Brune { 754ab8d36c9SPeter Brune SNES_FAS *fas; 7555fd66863SKarl Rupp 756ab8d36c9SPeter Brune PetscFunctionBegin; 757ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 758ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 759ab8d36c9SPeter Brune *mat = fas->inject; 760ab8d36c9SPeter Brune PetscFunctionReturn(0); 761ab8d36c9SPeter Brune } 762ab8d36c9SPeter Brune 763ab8d36c9SPeter Brune #undef __FUNCT__ 764ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRScale" 765ab8d36c9SPeter Brune /*@ 766ab8d36c9SPeter Brune SNESFASCycleGetRScale - Gets the injection on this level 767ab8d36c9SPeter Brune 768ab8d36c9SPeter Brune Logically Collective on SNES 769ab8d36c9SPeter Brune 770ab8d36c9SPeter Brune Input Parameters: 771ab8d36c9SPeter Brune . snes - the multigrid context 772ab8d36c9SPeter Brune 773ab8d36c9SPeter Brune Output Parameters: 774ab8d36c9SPeter Brune . mat - the restriction operator on this level 775ab8d36c9SPeter Brune 776ab8d36c9SPeter Brune Level: developer 777ab8d36c9SPeter Brune 778ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 779ab8d36c9SPeter Brune 780ab8d36c9SPeter Brune .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale() 781ab8d36c9SPeter Brune @*/ 782ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 783ab8d36c9SPeter Brune { 784ab8d36c9SPeter Brune SNES_FAS *fas; 7855fd66863SKarl Rupp 786ab8d36c9SPeter Brune PetscFunctionBegin; 787ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 788ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 789ab8d36c9SPeter Brune *vec = fas->rscale; 790ab8d36c9SPeter Brune PetscFunctionReturn(0); 791ab8d36c9SPeter Brune } 792ab8d36c9SPeter Brune 793ab8d36c9SPeter Brune #undef __FUNCT__ 794ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleIsFine" 795ab8d36c9SPeter Brune /*@ 796ab8d36c9SPeter Brune SNESFASCycleIsFine - Determines if a given cycle is the fine level. 797ab8d36c9SPeter Brune 798ab8d36c9SPeter Brune Logically Collective on SNES 799ab8d36c9SPeter Brune 800ab8d36c9SPeter Brune Input Parameters: 801ab8d36c9SPeter Brune . snes - the FAS context 802ab8d36c9SPeter Brune 803ab8d36c9SPeter Brune Output Parameters: 804ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not 805ab8d36c9SPeter Brune 806ab8d36c9SPeter Brune Level: advanced 807ab8d36c9SPeter Brune 808ab8d36c9SPeter Brune .keywords: SNES, FAS 809ab8d36c9SPeter Brune 810ab8d36c9SPeter Brune .seealso: SNESFASSetLevels() 811ab8d36c9SPeter Brune @*/ 812ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 813ab8d36c9SPeter Brune { 814ab8d36c9SPeter Brune SNES_FAS *fas; 8155fd66863SKarl Rupp 816ab8d36c9SPeter Brune PetscFunctionBegin; 817ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 818ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 8191aa26658SKarl Rupp if (fas->level == fas->levels - 1) *flg = PETSC_TRUE; 8201aa26658SKarl Rupp else *flg = PETSC_FALSE; 821ab8d36c9SPeter Brune PetscFunctionReturn(0); 822ab8d36c9SPeter Brune } 823ab8d36c9SPeter Brune 824ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */ 825ab8d36c9SPeter Brune 826ab8d36c9SPeter Brune #undef __FUNCT__ 827ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 828ab8d36c9SPeter Brune /*@ 829ab8d36c9SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 830ab8d36c9SPeter Brune interpolation from l-1 to the lth level 831ab8d36c9SPeter Brune 832ab8d36c9SPeter Brune Input Parameters: 833ab8d36c9SPeter Brune + snes - the multigrid context 834ab8d36c9SPeter Brune . mat - the interpolation operator 835ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 836ab8d36c9SPeter Brune 837ab8d36c9SPeter Brune Level: advanced 838ab8d36c9SPeter Brune 839ab8d36c9SPeter Brune Notes: 840ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction 841ab8d36c9SPeter Brune for the same level. 842ab8d36c9SPeter Brune 843ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 844ab8d36c9SPeter Brune out from the matrix size which one it is. 845ab8d36c9SPeter Brune 846ab8d36c9SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 847ab8d36c9SPeter Brune 848ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 849ab8d36c9SPeter Brune @*/ 8500adebc6cSBarry Smith PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) 8510adebc6cSBarry Smith { 85222d28d08SBarry Smith SNES_FAS *fas; 853ab8d36c9SPeter Brune PetscErrorCode ierr; 854ab8d36c9SPeter Brune SNES levelsnes; 85522d28d08SBarry Smith 856ab8d36c9SPeter Brune PetscFunctionBegin; 857ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 858ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 859ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 860ab8d36c9SPeter Brune ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 8611aa26658SKarl Rupp 862ab8d36c9SPeter Brune fas->interpolate = mat; 863ab8d36c9SPeter Brune PetscFunctionReturn(0); 864ab8d36c9SPeter Brune } 865ab8d36c9SPeter Brune 866ab8d36c9SPeter Brune #undef __FUNCT__ 867ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInterpolation" 868ab8d36c9SPeter Brune /*@ 869ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the 870ab8d36c9SPeter Brune interpolation from l-1 to the lth level 871ab8d36c9SPeter Brune 872ab8d36c9SPeter Brune Input Parameters: 873ab8d36c9SPeter Brune + snes - the multigrid context 874ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 875ab8d36c9SPeter Brune 876ab8d36c9SPeter Brune Output Parameters: 877ab8d36c9SPeter Brune . mat - the interpolation operator 878ab8d36c9SPeter Brune 879ab8d36c9SPeter Brune Level: advanced 880ab8d36c9SPeter Brune 881ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, interpolate, level 882ab8d36c9SPeter Brune 883ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale() 884ab8d36c9SPeter Brune @*/ 88522d28d08SBarry Smith PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) 88622d28d08SBarry Smith { 88722d28d08SBarry Smith SNES_FAS *fas; 888ab8d36c9SPeter Brune PetscErrorCode ierr; 889ab8d36c9SPeter Brune SNES levelsnes; 89022d28d08SBarry Smith 891ab8d36c9SPeter Brune PetscFunctionBegin; 892ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 893ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 894ab8d36c9SPeter Brune *mat = fas->interpolate; 895ab8d36c9SPeter Brune PetscFunctionReturn(0); 896ab8d36c9SPeter Brune } 897ab8d36c9SPeter Brune 898ab8d36c9SPeter Brune #undef __FUNCT__ 899ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 900ab8d36c9SPeter Brune /*@ 901ab8d36c9SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 902ab8d36c9SPeter Brune from level l to l-1. 903ab8d36c9SPeter Brune 904ab8d36c9SPeter Brune Input Parameters: 905ab8d36c9SPeter Brune + snes - the multigrid context 906ab8d36c9SPeter Brune . mat - the restriction matrix 907ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 908ab8d36c9SPeter Brune 909ab8d36c9SPeter Brune Level: advanced 910ab8d36c9SPeter Brune 911ab8d36c9SPeter Brune Notes: 912ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation 913ab8d36c9SPeter Brune for the same level. 914ab8d36c9SPeter Brune 915ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 916ab8d36c9SPeter Brune out from the matrix size which one it is. 917ab8d36c9SPeter Brune 918ab8d36c9SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 919ab8d36c9SPeter Brune is used. 920ab8d36c9SPeter Brune 921ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 922ab8d36c9SPeter Brune 923ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 924ab8d36c9SPeter Brune @*/ 92522d28d08SBarry Smith PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) 92622d28d08SBarry Smith { 92722d28d08SBarry Smith SNES_FAS *fas; 928ab8d36c9SPeter Brune PetscErrorCode ierr; 929ab8d36c9SPeter Brune SNES levelsnes; 93022d28d08SBarry Smith 931ab8d36c9SPeter Brune PetscFunctionBegin; 932ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 933ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 934ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 935ab8d36c9SPeter Brune ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 9361aa26658SKarl Rupp 937ab8d36c9SPeter Brune fas->restrct = mat; 938ab8d36c9SPeter Brune PetscFunctionReturn(0); 939ab8d36c9SPeter Brune } 940ab8d36c9SPeter Brune 941ab8d36c9SPeter Brune #undef __FUNCT__ 942ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetRestriction" 943ab8d36c9SPeter Brune /*@ 944ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the 945ab8d36c9SPeter Brune restriction from l to the l-1th level 946ab8d36c9SPeter Brune 947ab8d36c9SPeter Brune Input Parameters: 948ab8d36c9SPeter Brune + snes - the multigrid context 949ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 950ab8d36c9SPeter Brune 951ab8d36c9SPeter Brune Output Parameters: 952ab8d36c9SPeter Brune . mat - the interpolation operator 953ab8d36c9SPeter Brune 954ab8d36c9SPeter Brune Level: advanced 955ab8d36c9SPeter Brune 956ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, restrict, level 957ab8d36c9SPeter Brune 958ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale() 959ab8d36c9SPeter Brune @*/ 96022d28d08SBarry Smith PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) 96122d28d08SBarry Smith { 96222d28d08SBarry Smith SNES_FAS *fas; 963ab8d36c9SPeter Brune PetscErrorCode ierr; 964ab8d36c9SPeter Brune SNES levelsnes; 96522d28d08SBarry Smith 966ab8d36c9SPeter Brune PetscFunctionBegin; 967ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 968ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 969ab8d36c9SPeter Brune *mat = fas->restrct; 970ab8d36c9SPeter Brune PetscFunctionReturn(0); 971ab8d36c9SPeter Brune } 972ab8d36c9SPeter Brune 973ab8d36c9SPeter Brune 974ab8d36c9SPeter Brune #undef __FUNCT__ 975ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInjection" 976ab8d36c9SPeter Brune /*@ 977ab8d36c9SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 978ab8d36c9SPeter Brune from level l to l-1. 979ab8d36c9SPeter Brune 980ab8d36c9SPeter Brune Input Parameters: 981ab8d36c9SPeter Brune + snes - the multigrid context 982ab8d36c9SPeter Brune . mat - the restriction matrix 983ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 984ab8d36c9SPeter Brune 985ab8d36c9SPeter Brune Level: advanced 986ab8d36c9SPeter Brune 987ab8d36c9SPeter Brune Notes: 988ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to 989ab8d36c9SPeter Brune project the solution instead. 990ab8d36c9SPeter Brune 991ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 992ab8d36c9SPeter Brune 993ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 994ab8d36c9SPeter Brune @*/ 99522d28d08SBarry Smith PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) 99622d28d08SBarry Smith { 99722d28d08SBarry Smith SNES_FAS *fas; 998ab8d36c9SPeter Brune PetscErrorCode ierr; 999ab8d36c9SPeter Brune SNES levelsnes; 100022d28d08SBarry Smith 1001ab8d36c9SPeter Brune PetscFunctionBegin; 1002ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1003ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1004ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 1005ab8d36c9SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 10061aa26658SKarl Rupp 1007ab8d36c9SPeter Brune fas->inject = mat; 1008ab8d36c9SPeter Brune PetscFunctionReturn(0); 1009ab8d36c9SPeter Brune } 1010ab8d36c9SPeter Brune 1011ab8d36c9SPeter Brune 1012ab8d36c9SPeter Brune #undef __FUNCT__ 1013ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInjection" 1014ab8d36c9SPeter Brune /*@ 1015ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the 1016ab8d36c9SPeter Brune injection from l-1 to the lth level 1017ab8d36c9SPeter Brune 1018ab8d36c9SPeter Brune Input Parameters: 1019ab8d36c9SPeter Brune + snes - the multigrid context 1020ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 1021ab8d36c9SPeter Brune 1022ab8d36c9SPeter Brune Output Parameters: 1023ab8d36c9SPeter Brune . mat - the injection operator 1024ab8d36c9SPeter Brune 1025ab8d36c9SPeter Brune Level: advanced 1026ab8d36c9SPeter Brune 1027ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, restrict, level 1028ab8d36c9SPeter Brune 1029ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale() 1030ab8d36c9SPeter Brune @*/ 103122d28d08SBarry Smith PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) 103222d28d08SBarry Smith { 103322d28d08SBarry Smith SNES_FAS *fas; 1034ab8d36c9SPeter Brune PetscErrorCode ierr; 1035ab8d36c9SPeter Brune SNES levelsnes; 103622d28d08SBarry Smith 1037ab8d36c9SPeter Brune PetscFunctionBegin; 1038ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1039ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1040ab8d36c9SPeter Brune *mat = fas->inject; 1041ab8d36c9SPeter Brune PetscFunctionReturn(0); 1042ab8d36c9SPeter Brune } 1043ab8d36c9SPeter Brune 1044ab8d36c9SPeter Brune #undef __FUNCT__ 1045ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 1046ab8d36c9SPeter Brune /*@ 1047ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 1048ab8d36c9SPeter Brune operator from level l to l-1. 1049ab8d36c9SPeter Brune 1050ab8d36c9SPeter Brune Input Parameters: 1051ab8d36c9SPeter Brune + snes - the multigrid context 1052ab8d36c9SPeter Brune . rscale - the restriction scaling 1053ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 1054ab8d36c9SPeter Brune 1055ab8d36c9SPeter Brune Level: advanced 1056ab8d36c9SPeter Brune 1057ab8d36c9SPeter Brune Notes: 1058ab8d36c9SPeter Brune This is only used in the case that the injection is not set. 1059ab8d36c9SPeter Brune 1060ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 1061ab8d36c9SPeter Brune 1062ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1063ab8d36c9SPeter Brune @*/ 106422d28d08SBarry Smith PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) 106522d28d08SBarry Smith { 1066ab8d36c9SPeter Brune SNES_FAS *fas; 1067ab8d36c9SPeter Brune PetscErrorCode ierr; 1068ab8d36c9SPeter Brune SNES levelsnes; 106922d28d08SBarry Smith 1070ab8d36c9SPeter Brune PetscFunctionBegin; 1071ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1072ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1073ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 1074ab8d36c9SPeter Brune ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 10751aa26658SKarl Rupp 1076ab8d36c9SPeter Brune fas->rscale = rscale; 1077ab8d36c9SPeter Brune PetscFunctionReturn(0); 1078ab8d36c9SPeter Brune } 1079ab8d36c9SPeter Brune 1080ab8d36c9SPeter Brune #undef __FUNCT__ 1081ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmoother" 1082ab8d36c9SPeter Brune /*@ 1083ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level. 1084ab8d36c9SPeter Brune 1085ab8d36c9SPeter Brune Input Parameters: 1086ab8d36c9SPeter Brune + snes - the multigrid context 1087ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1088ab8d36c9SPeter Brune 1089ab8d36c9SPeter Brune Output Parameters: 1090ab8d36c9SPeter Brune smooth - the smoother 1091ab8d36c9SPeter Brune 1092ab8d36c9SPeter Brune Level: advanced 1093ab8d36c9SPeter Brune 1094ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1095ab8d36c9SPeter Brune 1096ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1097ab8d36c9SPeter Brune @*/ 109822d28d08SBarry Smith PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) 109922d28d08SBarry Smith { 1100ab8d36c9SPeter Brune SNES_FAS *fas; 1101ab8d36c9SPeter Brune PetscErrorCode ierr; 1102ab8d36c9SPeter Brune SNES levelsnes; 110322d28d08SBarry Smith 1104ab8d36c9SPeter Brune PetscFunctionBegin; 1105ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1106ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1107ab8d36c9SPeter Brune if (!fas->smoothd) { 1108ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 1109ab8d36c9SPeter Brune } 1110ab8d36c9SPeter Brune *smooth = fas->smoothd; 1111ab8d36c9SPeter Brune PetscFunctionReturn(0); 1112ab8d36c9SPeter Brune } 1113ab8d36c9SPeter Brune 1114ab8d36c9SPeter Brune #undef __FUNCT__ 1115ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherDown" 1116ab8d36c9SPeter Brune /*@ 1117ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level. 1118ab8d36c9SPeter Brune 1119ab8d36c9SPeter Brune Input Parameters: 1120ab8d36c9SPeter Brune + snes - the multigrid context 1121ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 1122ab8d36c9SPeter Brune 1123ab8d36c9SPeter Brune Output Parameters: 1124ab8d36c9SPeter Brune smooth - the smoother 1125ab8d36c9SPeter Brune 1126ab8d36c9SPeter Brune Level: advanced 1127ab8d36c9SPeter Brune 1128ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1129ab8d36c9SPeter Brune 1130ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1131ab8d36c9SPeter Brune @*/ 113222d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) 113322d28d08SBarry Smith { 1134ab8d36c9SPeter Brune SNES_FAS *fas; 1135ab8d36c9SPeter Brune PetscErrorCode ierr; 1136ab8d36c9SPeter Brune SNES levelsnes; 113722d28d08SBarry Smith 1138ab8d36c9SPeter Brune PetscFunctionBegin; 1139ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1140ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1141ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1142ab8d36c9SPeter Brune if (!fas->smoothd) { 1143ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 1144ab8d36c9SPeter Brune } 1145ab8d36c9SPeter Brune if (!fas->smoothu) { 1146ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 1147ab8d36c9SPeter Brune } 1148ab8d36c9SPeter Brune *smooth = fas->smoothd; 1149ab8d36c9SPeter Brune PetscFunctionReturn(0); 1150ab8d36c9SPeter Brune } 1151ab8d36c9SPeter Brune 1152ab8d36c9SPeter Brune #undef __FUNCT__ 1153ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherUp" 1154ab8d36c9SPeter Brune /*@ 1155ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level. 1156ab8d36c9SPeter Brune 1157ab8d36c9SPeter Brune Input Parameters: 1158ab8d36c9SPeter Brune + snes - the multigrid context 1159ab8d36c9SPeter Brune - level - the level (0 is coarsest) 1160ab8d36c9SPeter Brune 1161ab8d36c9SPeter Brune Output Parameters: 1162ab8d36c9SPeter Brune smooth - the smoother 1163ab8d36c9SPeter Brune 1164ab8d36c9SPeter Brune Level: advanced 1165ab8d36c9SPeter Brune 1166ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1167ab8d36c9SPeter Brune 1168ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1169ab8d36c9SPeter Brune @*/ 117022d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) 117122d28d08SBarry Smith { 1172ab8d36c9SPeter Brune SNES_FAS *fas; 1173ab8d36c9SPeter Brune PetscErrorCode ierr; 1174ab8d36c9SPeter Brune SNES levelsnes; 117522d28d08SBarry Smith 1176ab8d36c9SPeter Brune PetscFunctionBegin; 1177ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1178ab8d36c9SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1179ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1180ab8d36c9SPeter Brune if (!fas->smoothd) { 1181ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 1182ab8d36c9SPeter Brune } 1183ab8d36c9SPeter Brune if (!fas->smoothu) { 1184ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 1185ab8d36c9SPeter Brune } 1186ab8d36c9SPeter Brune *smooth = fas->smoothu; 1187ab8d36c9SPeter Brune PetscFunctionReturn(0); 1188ab8d36c9SPeter Brune } 1189d6ad1212SPeter Brune 1190d6ad1212SPeter Brune #undef __FUNCT__ 1191d6ad1212SPeter Brune #define __FUNCT__ "SNESFASGetCoarseSolve" 1192d6ad1212SPeter Brune /*@ 1193d6ad1212SPeter Brune SNESFASGetCoarseSolve - Gets the coarsest solver. 1194d6ad1212SPeter Brune 1195d6ad1212SPeter Brune Input Parameters: 1196d6ad1212SPeter Brune + snes - the multigrid context 1197d6ad1212SPeter Brune 1198d6ad1212SPeter Brune Output Parameters: 1199d6ad1212SPeter Brune solve - the coarse-level solver 1200d6ad1212SPeter Brune 1201d6ad1212SPeter Brune Level: advanced 1202d6ad1212SPeter Brune 1203d6ad1212SPeter Brune .keywords: FAS, MG, get, multigrid, solver, coarse 1204d6ad1212SPeter Brune 1205d6ad1212SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1206d6ad1212SPeter Brune @*/ 120722d28d08SBarry Smith PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *smooth) 120822d28d08SBarry Smith { 1209d6ad1212SPeter Brune SNES_FAS *fas; 1210d6ad1212SPeter Brune PetscErrorCode ierr; 1211d6ad1212SPeter Brune SNES levelsnes; 121222d28d08SBarry Smith 1213d6ad1212SPeter Brune PetscFunctionBegin; 1214d6ad1212SPeter Brune ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr); 1215d6ad1212SPeter Brune fas = (SNES_FAS*)levelsnes->data; 1216d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1217d6ad1212SPeter Brune if (!fas->smoothd) { 1218d6ad1212SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 1219d6ad1212SPeter Brune } 1220d6ad1212SPeter Brune *smooth = fas->smoothd; 1221d6ad1212SPeter Brune PetscFunctionReturn(0); 1222d6ad1212SPeter Brune } 1223928e959bSPeter Brune 1224928e959bSPeter Brune #undef __FUNCT__ 1225928e959bSPeter Brune #define __FUNCT__ "SNESFASFullSetDownSweep" 1226928e959bSPeter Brune /*@ 1227928e959bSPeter Brune SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS 1228928e959bSPeter Brune 1229928e959bSPeter Brune Logically Collective on SNES 1230928e959bSPeter Brune 1231928e959bSPeter Brune Input Parameters: 1232928e959bSPeter Brune + snes - the multigrid context 1233928e959bSPeter Brune - swp - whether to downsweep or not 1234928e959bSPeter Brune 1235928e959bSPeter Brune Options Database Key: 1236928e959bSPeter Brune . -snes_fas_full_downsweep - Sets number of pre-smoothing steps 1237928e959bSPeter Brune 1238928e959bSPeter Brune Level: advanced 1239928e959bSPeter Brune 1240928e959bSPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 1241928e959bSPeter Brune 1242928e959bSPeter Brune .seealso: SNESFASSetNumberSmoothUp() 1243928e959bSPeter Brune @*/ 1244928e959bSPeter Brune PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp) 1245928e959bSPeter Brune { 1246928e959bSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 1247928e959bSPeter Brune PetscErrorCode ierr = 0; 1248928e959bSPeter Brune 1249928e959bSPeter Brune PetscFunctionBegin; 1250928e959bSPeter Brune fas->full_downsweep = swp; 1251928e959bSPeter Brune if (fas->next) { 1252928e959bSPeter Brune ierr = SNESFASFullSetDownSweep(fas->next,swp);CHKERRQ(ierr); 1253928e959bSPeter Brune } 1254928e959bSPeter Brune PetscFunctionReturn(0); 1255928e959bSPeter Brune } 1256