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 17583a1250SSatish Balay - fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE 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; 27ab8d36c9SPeter Brune PetscFunctionBegin; 28ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 29ab8d36c9SPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 30ab8d36c9SPeter Brune fas->fastype = fastype; 31f7620de1SMatthew G Knepley if (fas->next) {ierr = SNESFASSetType(fas->next, fastype);CHKERRQ(ierr);} 32ab8d36c9SPeter Brune PetscFunctionReturn(0); 33ab8d36c9SPeter Brune } 34ab8d36c9SPeter Brune 35ab8d36c9SPeter Brune 36ab8d36c9SPeter Brune #undef __FUNCT__ 37ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetType" 38ab8d36c9SPeter Brune /*@ 39ab8d36c9SPeter Brune SNESFASGetType - Sets the update and correction type used for FAS. 40ab8d36c9SPeter Brune 41ab8d36c9SPeter Brune Logically Collective 42ab8d36c9SPeter Brune 43ab8d36c9SPeter Brune Input Parameters: 44ab8d36c9SPeter Brune . snes - FAS context 45ab8d36c9SPeter Brune 46ab8d36c9SPeter Brune Output Parameters: 47ab8d36c9SPeter Brune . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE 48ab8d36c9SPeter Brune 49583a1250SSatish Balay Level: intermediate 50583a1250SSatish Balay 51ab8d36c9SPeter Brune .seealso: PCMGSetType() 52ab8d36c9SPeter Brune @*/ 53ab8d36c9SPeter Brune PetscErrorCode SNESFASGetType(SNES snes,SNESFASType *fastype) 54ab8d36c9SPeter Brune { 55ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 56ab8d36c9SPeter Brune 57ab8d36c9SPeter Brune PetscFunctionBegin; 58ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 59ab8d36c9SPeter Brune PetscValidPointer(fastype, 2); 60ab8d36c9SPeter Brune *fastype = fas->fastype; 61ab8d36c9SPeter Brune PetscFunctionReturn(0); 62ab8d36c9SPeter Brune } 63ab8d36c9SPeter Brune 64ab8d36c9SPeter Brune #undef __FUNCT__ 65ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 66ab8d36c9SPeter Brune /*@C 67ab8d36c9SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 68ab8d36c9SPeter Brune Must be called before any other FAS routine. 69ab8d36c9SPeter Brune 70ab8d36c9SPeter Brune Input Parameters: 71ab8d36c9SPeter Brune + snes - the snes context 72ab8d36c9SPeter Brune . levels - the number of levels 73ab8d36c9SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 74ab8d36c9SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 75ab8d36c9SPeter Brune Fortran. 76ab8d36c9SPeter Brune 77ab8d36c9SPeter Brune Level: intermediate 78ab8d36c9SPeter Brune 79ab8d36c9SPeter Brune Notes: 80ab8d36c9SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 81ab8d36c9SPeter Brune for setting the level options rather than the -fas_coarse prefix. 82ab8d36c9SPeter Brune 83ab8d36c9SPeter Brune .keywords: FAS, MG, set, levels, multigrid 84ab8d36c9SPeter Brune 85ab8d36c9SPeter Brune .seealso: SNESFASGetLevels() 86ab8d36c9SPeter Brune @*/ 87ab8d36c9SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 88ab8d36c9SPeter Brune PetscErrorCode ierr; 89ab8d36c9SPeter Brune PetscInt i; 90ab8d36c9SPeter Brune const char *optionsprefix; 91ab8d36c9SPeter Brune char tprefix[128]; 92ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 93ab8d36c9SPeter Brune SNES prevsnes; 94ab8d36c9SPeter Brune MPI_Comm comm; 95ab8d36c9SPeter Brune PetscFunctionBegin; 96ab8d36c9SPeter Brune comm = ((PetscObject)snes)->comm; 97ab8d36c9SPeter Brune if (levels == fas->levels) { 98ab8d36c9SPeter Brune if (!comms) 99ab8d36c9SPeter Brune PetscFunctionReturn(0); 100ab8d36c9SPeter Brune } 101ab8d36c9SPeter Brune /* user has changed the number of levels; reset */ 102ab8d36c9SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 103ab8d36c9SPeter Brune /* destroy any coarser levels if necessary */ 104ab8d36c9SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 105ab8d36c9SPeter Brune fas->next = PETSC_NULL; 106ab8d36c9SPeter Brune fas->previous = PETSC_NULL; 107ab8d36c9SPeter Brune prevsnes = snes; 108ab8d36c9SPeter Brune /* setup the finest level */ 109ab8d36c9SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 110ab8d36c9SPeter Brune for (i = levels - 1; i >= 0; i--) { 111ab8d36c9SPeter Brune if (comms) comm = comms[i]; 112ab8d36c9SPeter Brune fas->level = i; 113ab8d36c9SPeter Brune fas->levels = levels; 114ab8d36c9SPeter Brune fas->fine = snes; 115ab8d36c9SPeter Brune fas->next = PETSC_NULL; 116ab8d36c9SPeter Brune if (i > 0) { 117ab8d36c9SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 118ab8d36c9SPeter Brune sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level); 119ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(fas->next,tprefix);CHKERRQ(ierr); 120ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(fas->next,optionsprefix);CHKERRQ(ierr); 121ab8d36c9SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 122ab8d36c9SPeter Brune ierr = SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);CHKERRQ(ierr); 123ab8d36c9SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 124ab8d36c9SPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 125ab8d36c9SPeter Brune prevsnes = fas->next; 126ab8d36c9SPeter Brune fas = (SNES_FAS *)prevsnes->data; 127ab8d36c9SPeter Brune } 128ab8d36c9SPeter Brune } 129ab8d36c9SPeter Brune PetscFunctionReturn(0); 130ab8d36c9SPeter Brune } 131ab8d36c9SPeter Brune 132ab8d36c9SPeter Brune 133ab8d36c9SPeter Brune #undef __FUNCT__ 134ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 135ab8d36c9SPeter Brune /*@ 136ab8d36c9SPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 137ab8d36c9SPeter Brune 138ab8d36c9SPeter Brune Input Parameter: 139ab8d36c9SPeter Brune . snes - the nonlinear solver context 140ab8d36c9SPeter Brune 141ab8d36c9SPeter Brune Output parameter: 142ab8d36c9SPeter Brune . levels - the number of levels 143ab8d36c9SPeter Brune 144ab8d36c9SPeter Brune Level: advanced 145ab8d36c9SPeter Brune 146ab8d36c9SPeter Brune .keywords: MG, get, levels, multigrid 147ab8d36c9SPeter Brune 148ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 149ab8d36c9SPeter Brune @*/ 150ab8d36c9SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 151ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 152ab8d36c9SPeter Brune PetscFunctionBegin; 153ab8d36c9SPeter Brune *levels = fas->levels; 154ab8d36c9SPeter Brune PetscFunctionReturn(0); 155ab8d36c9SPeter Brune } 156ab8d36c9SPeter Brune 157ab8d36c9SPeter Brune 158ab8d36c9SPeter Brune #undef __FUNCT__ 159ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetCycleSNES" 160ab8d36c9SPeter Brune /*@ 161ab8d36c9SPeter Brune SNESFASGetCycleSNES - Gets the SNES corresponding to a particular 162ab8d36c9SPeter Brune level of the FAS hierarchy. 163ab8d36c9SPeter Brune 164ab8d36c9SPeter Brune Input Parameters: 165ab8d36c9SPeter Brune + snes - the multigrid context 166ab8d36c9SPeter Brune level - the level to get 167ab8d36c9SPeter Brune - lsnes - whether to use the nonlinear smoother or not 168ab8d36c9SPeter Brune 169ab8d36c9SPeter Brune Level: advanced 170ab8d36c9SPeter Brune 171ab8d36c9SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 172ab8d36c9SPeter Brune 173ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 174ab8d36c9SPeter Brune @*/ 175ab8d36c9SPeter Brune PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes) { 176ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 177ab8d36c9SPeter Brune PetscInt i; 178ab8d36c9SPeter Brune 179ab8d36c9SPeter Brune PetscFunctionBegin; 180ab8d36c9SPeter Brune if (level > fas->levels-1) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels); 181ab8d36c9SPeter Brune if (fas->level != fas->levels - 1) 182ab8d36c9SPeter Brune SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level); 183ab8d36c9SPeter Brune 184ab8d36c9SPeter Brune *lsnes = snes; 185ab8d36c9SPeter Brune for (i = fas->level; i > level; i--) { 186ab8d36c9SPeter Brune *lsnes = fas->next; 187ab8d36c9SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 188ab8d36c9SPeter Brune } 189ab8d36c9SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 190ab8d36c9SPeter Brune PetscFunctionReturn(0); 191ab8d36c9SPeter Brune } 192ab8d36c9SPeter Brune 193ab8d36c9SPeter Brune #undef __FUNCT__ 194ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 195ab8d36c9SPeter Brune /*@ 196ab8d36c9SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 197ab8d36c9SPeter Brune use on all levels. 198ab8d36c9SPeter Brune 199ab8d36c9SPeter Brune Logically Collective on SNES 200ab8d36c9SPeter Brune 201ab8d36c9SPeter Brune Input Parameters: 202ab8d36c9SPeter Brune + snes - the multigrid context 203ab8d36c9SPeter Brune - n - the number of smoothing steps 204ab8d36c9SPeter Brune 205ab8d36c9SPeter Brune Options Database Key: 206ab8d36c9SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 207ab8d36c9SPeter Brune 208ab8d36c9SPeter Brune Level: advanced 209ab8d36c9SPeter Brune 210ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 211ab8d36c9SPeter Brune 212ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 213ab8d36c9SPeter Brune @*/ 214ab8d36c9SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 215ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 216ab8d36c9SPeter Brune PetscErrorCode ierr = 0; 217ab8d36c9SPeter Brune PetscFunctionBegin; 218ab8d36c9SPeter Brune fas->max_up_it = n; 219*656ede7eSPeter Brune if (!fas->smoothu && fas->level != 0) { 220ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu); 221ab8d36c9SPeter Brune } 222*656ede7eSPeter Brune if (fas->smoothu)ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr); 223ab8d36c9SPeter Brune if (fas->next) { 224ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 225ab8d36c9SPeter Brune } 226ab8d36c9SPeter Brune PetscFunctionReturn(0); 227ab8d36c9SPeter Brune } 228ab8d36c9SPeter Brune 229ab8d36c9SPeter Brune #undef __FUNCT__ 230ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 231ab8d36c9SPeter Brune /*@ 232ab8d36c9SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 233ab8d36c9SPeter Brune use on all levels. 234ab8d36c9SPeter Brune 235ab8d36c9SPeter Brune Logically Collective on SNES 236ab8d36c9SPeter Brune 237ab8d36c9SPeter Brune Input Parameters: 238ab8d36c9SPeter Brune + snes - the multigrid context 239ab8d36c9SPeter Brune - n - the number of smoothing steps 240ab8d36c9SPeter Brune 241ab8d36c9SPeter Brune Options Database Key: 242ab8d36c9SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 243ab8d36c9SPeter Brune 244ab8d36c9SPeter Brune Level: advanced 245ab8d36c9SPeter Brune 246ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 247ab8d36c9SPeter Brune 248ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 249ab8d36c9SPeter Brune @*/ 250ab8d36c9SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 251ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 252ab8d36c9SPeter Brune PetscErrorCode ierr = 0; 253ab8d36c9SPeter Brune PetscFunctionBegin; 254ab8d36c9SPeter Brune if (!fas->smoothu) { 255ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu); 256ab8d36c9SPeter Brune } 257ab8d36c9SPeter Brune if (!fas->smoothd) { 258ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd); 259ab8d36c9SPeter Brune } 260ab8d36c9SPeter Brune ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr); 261ab8d36c9SPeter Brune fas->max_down_it = n; 262ab8d36c9SPeter Brune if (fas->next) { 263ab8d36c9SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 264ab8d36c9SPeter Brune } 265ab8d36c9SPeter Brune PetscFunctionReturn(0); 266ab8d36c9SPeter Brune } 267ab8d36c9SPeter Brune 268ab8d36c9SPeter Brune 269ab8d36c9SPeter Brune #undef __FUNCT__ 270ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetCycles" 271ab8d36c9SPeter Brune /*@ 272ab8d36c9SPeter Brune SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 273ab8d36c9SPeter Brune complicated cycling. 274ab8d36c9SPeter Brune 275ab8d36c9SPeter Brune Logically Collective on SNES 276ab8d36c9SPeter Brune 277ab8d36c9SPeter Brune Input Parameters: 278ab8d36c9SPeter Brune + snes - the multigrid context 279ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 280ab8d36c9SPeter Brune 281ab8d36c9SPeter Brune Options Database Key: 282ab8d36c9SPeter Brune $ -snes_fas_cycles 1 or 2 283ab8d36c9SPeter Brune 284ab8d36c9SPeter Brune Level: advanced 285ab8d36c9SPeter Brune 286ab8d36c9SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 287ab8d36c9SPeter Brune 288ab8d36c9SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 289ab8d36c9SPeter Brune @*/ 290ab8d36c9SPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 291ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 292ab8d36c9SPeter Brune PetscErrorCode ierr; 293ab8d36c9SPeter Brune PetscBool isFine; 294ab8d36c9SPeter Brune PetscFunctionBegin; 295ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine); 296ab8d36c9SPeter Brune fas->n_cycles = cycles; 297ab8d36c9SPeter Brune if (!isFine) 298ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 299ab8d36c9SPeter Brune if (fas->next) { 300ab8d36c9SPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 301ab8d36c9SPeter Brune } 302ab8d36c9SPeter Brune PetscFunctionReturn(0); 303ab8d36c9SPeter Brune } 304ab8d36c9SPeter Brune 305c8c899caSPeter Brune 306c8c899caSPeter Brune #undef __FUNCT__ 307c8c899caSPeter Brune #define __FUNCT__ "SNESFASSetMonitor" 308c8c899caSPeter Brune /*@ 309c8c899caSPeter Brune SNESFASSetMonitor - Sets the method-specific cycle monitoring 310c8c899caSPeter Brune 311c8c899caSPeter Brune Logically Collective on SNES 312c8c899caSPeter Brune 313c8c899caSPeter Brune Input Parameters: 314c8c899caSPeter Brune + snes - the FAS context 315c8c899caSPeter Brune - flg - monitor or not 316c8c899caSPeter Brune 317c8c899caSPeter Brune Level: advanced 318c8c899caSPeter Brune 319c8c899caSPeter Brune .keywords: FAS, monitor 320c8c899caSPeter Brune 321c8c899caSPeter Brune .seealso: SNESFASSetCyclesOnLevel() 322c8c899caSPeter Brune @*/ 323c8c899caSPeter Brune PetscErrorCode SNESFASSetMonitor(SNES snes, PetscBool flg) { 324c8c899caSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 325c8c899caSPeter Brune PetscErrorCode ierr; 326c8c899caSPeter Brune PetscBool isFine; 327c8c899caSPeter Brune PetscInt i, levels = fas->levels; 328c8c899caSPeter Brune SNES levelsnes; 329c8c899caSPeter Brune PetscFunctionBegin; 330c8c899caSPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 331c8c899caSPeter Brune if (isFine) { 332c8c899caSPeter Brune for (i = 0; i < levels; i++) { 333c8c899caSPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes); 334c8c899caSPeter Brune fas = (SNES_FAS *)levelsnes->data; 335c8c899caSPeter Brune if (flg) { 336c8c899caSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)levelsnes)->comm);CHKERRQ(ierr); 337c8c899caSPeter Brune /* set the monitors for the upsmoother and downsmoother */ 338c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 339c8c899caSPeter Brune ierr = SNESMonitorSet(levelsnes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 340c8c899caSPeter Brune } else { 341c8c899caSPeter Brune /* unset the monitors on the coarse levels */ 342c8c899caSPeter Brune if (i != fas->levels - 1) { 343c8c899caSPeter Brune ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr); 344c8c899caSPeter Brune } 345c8c899caSPeter Brune } 346c8c899caSPeter Brune } 347c8c899caSPeter Brune } 348c8c899caSPeter Brune PetscFunctionReturn(0); 349c8c899caSPeter Brune } 350c8c899caSPeter Brune 351ab8d36c9SPeter Brune #undef __FUNCT__ 352ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleCreateSmoother_Private" 353ab8d36c9SPeter Brune /* 354ab8d36c9SPeter Brune Creates the default smoother type. 355ab8d36c9SPeter Brune 356ab8d36c9SPeter Brune This is SNESNRICHARDSON on each fine level and SNESLS on the coarse level. 357ab8d36c9SPeter Brune 358ab8d36c9SPeter Brune */ 359ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) { 360ab8d36c9SPeter Brune SNES_FAS *fas; 361ab8d36c9SPeter Brune const char *optionsprefix; 362ab8d36c9SPeter Brune char tprefix[128]; 363ab8d36c9SPeter Brune PetscErrorCode ierr; 364ab8d36c9SPeter Brune SNES nsmooth; 365ab8d36c9SPeter Brune PetscFunctionBegin; 366ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 367ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 368ab8d36c9SPeter Brune ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr); 369ab8d36c9SPeter Brune /* create the default smoother */ 370ab8d36c9SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &nsmooth);CHKERRQ(ierr); 371ab8d36c9SPeter Brune if (fas->level == 0) { 372ab8d36c9SPeter Brune sprintf(tprefix,"fas_coarse_"); 373ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 374ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 375ab8d36c9SPeter Brune ierr = SNESSetType(nsmooth, SNESLS);CHKERRQ(ierr); 376e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr); 377ab8d36c9SPeter Brune } else { 378ab8d36c9SPeter Brune sprintf(tprefix,"fas_levels_%d_",(int)fas->level); 379ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr); 380ab8d36c9SPeter Brune ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr); 381ab8d36c9SPeter Brune ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr); 382e70c42e5SPeter Brune ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr); 383ab8d36c9SPeter Brune } 384ab8d36c9SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 385f89ba88eSPeter Brune ierr = PetscLogObjectParent(snes,nsmooth);CHKERRQ(ierr); 386f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr); 387ab8d36c9SPeter Brune *smooth = nsmooth; 388ab8d36c9SPeter Brune PetscFunctionReturn(0); 389ab8d36c9SPeter Brune } 390ab8d36c9SPeter Brune 391ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */ 392ab8d36c9SPeter Brune 393ab8d36c9SPeter Brune #undef __FUNCT__ 394ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleSetCycles" 395ab8d36c9SPeter Brune /*@ 396ab8d36c9SPeter Brune SNESFASCycleSetCycles - Sets the number of cycles on a particular level. 397ab8d36c9SPeter Brune 398ab8d36c9SPeter Brune Logically Collective on SNES 399ab8d36c9SPeter Brune 400ab8d36c9SPeter Brune Input Parameters: 401ab8d36c9SPeter Brune + snes - the multigrid context 402ab8d36c9SPeter Brune . level - the level to set the number of cycles on 403ab8d36c9SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 404ab8d36c9SPeter Brune 405ab8d36c9SPeter Brune Level: advanced 406ab8d36c9SPeter Brune 407ab8d36c9SPeter Brune .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid 408ab8d36c9SPeter Brune 409ab8d36c9SPeter Brune .seealso: SNESFASSetCycles() 410ab8d36c9SPeter Brune @*/ 411ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) { 412ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 413ab8d36c9SPeter Brune PetscErrorCode ierr; 414ab8d36c9SPeter Brune PetscFunctionBegin; 415ab8d36c9SPeter Brune fas->n_cycles = cycles; 416ab8d36c9SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr); 417ab8d36c9SPeter Brune PetscFunctionReturn(0); 418ab8d36c9SPeter Brune } 419ab8d36c9SPeter Brune 420ab8d36c9SPeter Brune 421ab8d36c9SPeter Brune #undef __FUNCT__ 422ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmoother" 423ab8d36c9SPeter Brune /*@ 424ab8d36c9SPeter Brune SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level. 425ab8d36c9SPeter Brune 426ab8d36c9SPeter Brune Logically Collective on SNES 427ab8d36c9SPeter Brune 428ab8d36c9SPeter Brune Input Parameters: 429ab8d36c9SPeter Brune . snes - the multigrid context 430ab8d36c9SPeter Brune 431ab8d36c9SPeter Brune Output Parameters: 432ab8d36c9SPeter Brune . smooth - the smoother 433ab8d36c9SPeter Brune 434ab8d36c9SPeter Brune Level: advanced 435ab8d36c9SPeter Brune 436ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 437ab8d36c9SPeter Brune 438ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown() 439ab8d36c9SPeter Brune @*/ 440ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) 441ab8d36c9SPeter Brune { 442ab8d36c9SPeter Brune SNES_FAS *fas; 443ab8d36c9SPeter Brune PetscFunctionBegin; 444ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 445ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 446ab8d36c9SPeter Brune *smooth = fas->smoothd; 447ab8d36c9SPeter Brune PetscFunctionReturn(0); 448ab8d36c9SPeter Brune } 449ab8d36c9SPeter Brune #undef __FUNCT__ 450ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherUp" 451ab8d36c9SPeter Brune /*@ 452ab8d36c9SPeter Brune SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level. 453ab8d36c9SPeter Brune 454ab8d36c9SPeter Brune Logically Collective on SNES 455ab8d36c9SPeter Brune 456ab8d36c9SPeter Brune Input Parameters: 457ab8d36c9SPeter Brune . snes - the multigrid context 458ab8d36c9SPeter Brune 459ab8d36c9SPeter Brune Output Parameters: 460ab8d36c9SPeter Brune . smoothu - the smoother 461ab8d36c9SPeter Brune 462ab8d36c9SPeter Brune Notes: 463ab8d36c9SPeter Brune Returns the downsmoother if no up smoother is available. This enables transparent 464ab8d36c9SPeter Brune default behavior in the process of the solve. 465ab8d36c9SPeter Brune 466ab8d36c9SPeter Brune Level: advanced 467ab8d36c9SPeter Brune 468ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 469ab8d36c9SPeter Brune 470ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown() 471ab8d36c9SPeter Brune @*/ 472ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) 473ab8d36c9SPeter Brune { 474ab8d36c9SPeter Brune SNES_FAS *fas; 475ab8d36c9SPeter Brune PetscFunctionBegin; 476ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 477ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 478ab8d36c9SPeter Brune if (!fas->smoothu) { 479ab8d36c9SPeter Brune *smoothu = fas->smoothd; 480ab8d36c9SPeter Brune } else { 481ab8d36c9SPeter Brune *smoothu = fas->smoothu; 482ab8d36c9SPeter Brune } 483ab8d36c9SPeter Brune PetscFunctionReturn(0); 484ab8d36c9SPeter Brune } 485ab8d36c9SPeter Brune 486ab8d36c9SPeter Brune #undef __FUNCT__ 487ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherDown" 488ab8d36c9SPeter Brune /*@ 489ab8d36c9SPeter Brune SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level. 490ab8d36c9SPeter Brune 491ab8d36c9SPeter Brune Logically Collective on SNES 492ab8d36c9SPeter Brune 493ab8d36c9SPeter Brune Input Parameters: 494ab8d36c9SPeter Brune . snes - the multigrid context 495ab8d36c9SPeter Brune 496ab8d36c9SPeter Brune Output Parameters: 497ab8d36c9SPeter Brune . smoothd - the smoother 498ab8d36c9SPeter Brune 499ab8d36c9SPeter Brune Level: advanced 500ab8d36c9SPeter Brune 501ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 502ab8d36c9SPeter Brune 503ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 504ab8d36c9SPeter Brune @*/ 505ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) 506ab8d36c9SPeter Brune { 507ab8d36c9SPeter Brune SNES_FAS *fas; 508ab8d36c9SPeter Brune PetscFunctionBegin; 509ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 510ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 511ab8d36c9SPeter Brune *smoothd = fas->smoothd; 512ab8d36c9SPeter Brune PetscFunctionReturn(0); 513ab8d36c9SPeter Brune } 514ab8d36c9SPeter Brune 515ab8d36c9SPeter Brune 516ab8d36c9SPeter Brune #undef __FUNCT__ 517ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetCorrection" 518ab8d36c9SPeter Brune /*@ 519ab8d36c9SPeter Brune SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level 520ab8d36c9SPeter Brune 521ab8d36c9SPeter Brune Logically Collective on SNES 522ab8d36c9SPeter Brune 523ab8d36c9SPeter Brune Input Parameters: 524ab8d36c9SPeter Brune . snes - the multigrid context 525ab8d36c9SPeter Brune 526ab8d36c9SPeter Brune Output Parameters: 527ab8d36c9SPeter Brune . correction - the coarse correction on this level 528ab8d36c9SPeter Brune 529ab8d36c9SPeter Brune Notes: 530ab8d36c9SPeter Brune Returns PETSC_NULL on the coarsest level. 531ab8d36c9SPeter Brune 532ab8d36c9SPeter Brune Level: advanced 533ab8d36c9SPeter Brune 534ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 535ab8d36c9SPeter Brune 536ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 537ab8d36c9SPeter Brune @*/ 538ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) 539ab8d36c9SPeter Brune { 540ab8d36c9SPeter Brune SNES_FAS *fas; 541ab8d36c9SPeter Brune PetscFunctionBegin; 542ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 543ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 544ab8d36c9SPeter Brune *correction = fas->next; 545ab8d36c9SPeter Brune PetscFunctionReturn(0); 546ab8d36c9SPeter Brune } 547ab8d36c9SPeter Brune 548ab8d36c9SPeter Brune #undef __FUNCT__ 549ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInterpolation" 550ab8d36c9SPeter Brune /*@ 551ab8d36c9SPeter Brune SNESFASCycleGetInterpolation - Gets the interpolation on this level 552ab8d36c9SPeter Brune 553ab8d36c9SPeter Brune Logically Collective on SNES 554ab8d36c9SPeter Brune 555ab8d36c9SPeter Brune Input Parameters: 556ab8d36c9SPeter Brune . snes - the multigrid context 557ab8d36c9SPeter Brune 558ab8d36c9SPeter Brune Output Parameters: 559ab8d36c9SPeter Brune . mat - the interpolation operator on this level 560ab8d36c9SPeter Brune 561ab8d36c9SPeter Brune Level: developer 562ab8d36c9SPeter Brune 563ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 564ab8d36c9SPeter Brune 565ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother() 566ab8d36c9SPeter Brune @*/ 567ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) 568ab8d36c9SPeter Brune { 569ab8d36c9SPeter Brune SNES_FAS *fas; 570ab8d36c9SPeter Brune PetscFunctionBegin; 571ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 572ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 573ab8d36c9SPeter Brune *mat = fas->interpolate; 574ab8d36c9SPeter Brune PetscFunctionReturn(0); 575ab8d36c9SPeter Brune } 576ab8d36c9SPeter Brune 577ab8d36c9SPeter Brune 578ab8d36c9SPeter Brune #undef __FUNCT__ 579ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRestriction" 580ab8d36c9SPeter Brune /*@ 581ab8d36c9SPeter Brune SNESFASCycleGetRestriction - Gets the restriction on this level 582ab8d36c9SPeter Brune 583ab8d36c9SPeter Brune Logically Collective on SNES 584ab8d36c9SPeter Brune 585ab8d36c9SPeter Brune Input Parameters: 586ab8d36c9SPeter Brune . snes - the multigrid context 587ab8d36c9SPeter Brune 588ab8d36c9SPeter Brune Output Parameters: 589ab8d36c9SPeter Brune . mat - the restriction operator on this level 590ab8d36c9SPeter Brune 591ab8d36c9SPeter Brune Level: developer 592ab8d36c9SPeter Brune 593ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 594ab8d36c9SPeter Brune 595ab8d36c9SPeter Brune .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation() 596ab8d36c9SPeter Brune @*/ 597ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) 598ab8d36c9SPeter Brune { 599ab8d36c9SPeter Brune SNES_FAS *fas; 600ab8d36c9SPeter Brune PetscFunctionBegin; 601ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 602ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 603ab8d36c9SPeter Brune *mat = fas->restrct; 604ab8d36c9SPeter Brune PetscFunctionReturn(0); 605ab8d36c9SPeter Brune } 606ab8d36c9SPeter Brune 607ab8d36c9SPeter Brune 608ab8d36c9SPeter Brune #undef __FUNCT__ 609ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInjection" 610ab8d36c9SPeter Brune /*@ 611ab8d36c9SPeter Brune SNESFASCycleGetInjection - Gets the injection on this level 612ab8d36c9SPeter Brune 613ab8d36c9SPeter Brune Logically Collective on SNES 614ab8d36c9SPeter Brune 615ab8d36c9SPeter Brune Input Parameters: 616ab8d36c9SPeter Brune . snes - the multigrid context 617ab8d36c9SPeter Brune 618ab8d36c9SPeter Brune Output Parameters: 619ab8d36c9SPeter Brune . mat - the restriction operator on this level 620ab8d36c9SPeter Brune 621ab8d36c9SPeter Brune Level: developer 622ab8d36c9SPeter Brune 623ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 624ab8d36c9SPeter Brune 625ab8d36c9SPeter Brune .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction() 626ab8d36c9SPeter Brune @*/ 627ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) 628ab8d36c9SPeter Brune { 629ab8d36c9SPeter Brune SNES_FAS *fas; 630ab8d36c9SPeter Brune PetscFunctionBegin; 631ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 632ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 633ab8d36c9SPeter Brune *mat = fas->inject; 634ab8d36c9SPeter Brune PetscFunctionReturn(0); 635ab8d36c9SPeter Brune } 636ab8d36c9SPeter Brune 637ab8d36c9SPeter Brune #undef __FUNCT__ 638ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRScale" 639ab8d36c9SPeter Brune /*@ 640ab8d36c9SPeter Brune SNESFASCycleGetRScale - Gets the injection on this level 641ab8d36c9SPeter Brune 642ab8d36c9SPeter Brune Logically Collective on SNES 643ab8d36c9SPeter Brune 644ab8d36c9SPeter Brune Input Parameters: 645ab8d36c9SPeter Brune . snes - the multigrid context 646ab8d36c9SPeter Brune 647ab8d36c9SPeter Brune Output Parameters: 648ab8d36c9SPeter Brune . mat - the restriction operator on this level 649ab8d36c9SPeter Brune 650ab8d36c9SPeter Brune Level: developer 651ab8d36c9SPeter Brune 652ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid 653ab8d36c9SPeter Brune 654ab8d36c9SPeter Brune .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale() 655ab8d36c9SPeter Brune @*/ 656ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) 657ab8d36c9SPeter Brune { 658ab8d36c9SPeter Brune SNES_FAS *fas; 659ab8d36c9SPeter Brune PetscFunctionBegin; 660ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 661ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 662ab8d36c9SPeter Brune *vec = fas->rscale; 663ab8d36c9SPeter Brune PetscFunctionReturn(0); 664ab8d36c9SPeter Brune } 665ab8d36c9SPeter Brune 666ab8d36c9SPeter Brune #undef __FUNCT__ 667ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleIsFine" 668ab8d36c9SPeter Brune /*@ 669ab8d36c9SPeter Brune SNESFASCycleIsFine - Determines if a given cycle is the fine level. 670ab8d36c9SPeter Brune 671ab8d36c9SPeter Brune Logically Collective on SNES 672ab8d36c9SPeter Brune 673ab8d36c9SPeter Brune Input Parameters: 674ab8d36c9SPeter Brune . snes - the FAS context 675ab8d36c9SPeter Brune 676ab8d36c9SPeter Brune Output Parameters: 677ab8d36c9SPeter Brune . flg - indicates if this is the fine level or not 678ab8d36c9SPeter Brune 679ab8d36c9SPeter Brune Level: advanced 680ab8d36c9SPeter Brune 681ab8d36c9SPeter Brune .keywords: SNES, FAS 682ab8d36c9SPeter Brune 683ab8d36c9SPeter Brune .seealso: SNESFASSetLevels() 684ab8d36c9SPeter Brune @*/ 685ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) 686ab8d36c9SPeter Brune { 687ab8d36c9SPeter Brune SNES_FAS *fas; 688ab8d36c9SPeter Brune PetscFunctionBegin; 689ab8d36c9SPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 690ab8d36c9SPeter Brune fas = (SNES_FAS*)snes->data; 691ab8d36c9SPeter Brune if (fas->level == fas->levels - 1) { 692ab8d36c9SPeter Brune *flg = PETSC_TRUE; 693ab8d36c9SPeter Brune } else { 694ab8d36c9SPeter Brune *flg = PETSC_FALSE; 695ab8d36c9SPeter Brune } 696ab8d36c9SPeter Brune PetscFunctionReturn(0); 697ab8d36c9SPeter Brune } 698ab8d36c9SPeter Brune 699ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */ 700ab8d36c9SPeter Brune 701ab8d36c9SPeter Brune #undef __FUNCT__ 702ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 703ab8d36c9SPeter Brune /*@ 704ab8d36c9SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 705ab8d36c9SPeter Brune interpolation from l-1 to the lth level 706ab8d36c9SPeter Brune 707ab8d36c9SPeter Brune Input Parameters: 708ab8d36c9SPeter Brune + snes - the multigrid context 709ab8d36c9SPeter Brune . mat - the interpolation operator 710ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 711ab8d36c9SPeter Brune 712ab8d36c9SPeter Brune Level: advanced 713ab8d36c9SPeter Brune 714ab8d36c9SPeter Brune Notes: 715ab8d36c9SPeter Brune Usually this is the same matrix used also to set the restriction 716ab8d36c9SPeter Brune for the same level. 717ab8d36c9SPeter Brune 718ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 719ab8d36c9SPeter Brune out from the matrix size which one it is. 720ab8d36c9SPeter Brune 721ab8d36c9SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 722ab8d36c9SPeter Brune 723ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 724ab8d36c9SPeter Brune @*/ 725ab8d36c9SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 726ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 727ab8d36c9SPeter Brune PetscErrorCode ierr; 728ab8d36c9SPeter Brune SNES levelsnes; 729ab8d36c9SPeter Brune PetscFunctionBegin; 730ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 731ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 732ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 733ab8d36c9SPeter Brune ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 734ab8d36c9SPeter Brune fas->interpolate = mat; 735ab8d36c9SPeter Brune PetscFunctionReturn(0); 736ab8d36c9SPeter Brune } 737ab8d36c9SPeter Brune 738ab8d36c9SPeter Brune #undef __FUNCT__ 739ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInterpolation" 740ab8d36c9SPeter Brune /*@ 741ab8d36c9SPeter Brune SNESFASGetInterpolation - Gets the matrix used to calculate the 742ab8d36c9SPeter Brune interpolation from l-1 to the lth level 743ab8d36c9SPeter Brune 744ab8d36c9SPeter Brune Input Parameters: 745ab8d36c9SPeter Brune + snes - the multigrid context 746ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 747ab8d36c9SPeter Brune 748ab8d36c9SPeter Brune Output Parameters: 749ab8d36c9SPeter Brune . mat - the interpolation operator 750ab8d36c9SPeter Brune 751ab8d36c9SPeter Brune Level: advanced 752ab8d36c9SPeter Brune 753ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, interpolate, level 754ab8d36c9SPeter Brune 755ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale() 756ab8d36c9SPeter Brune @*/ 757ab8d36c9SPeter Brune PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) { 758ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 759ab8d36c9SPeter Brune PetscErrorCode ierr; 760ab8d36c9SPeter Brune SNES levelsnes; 761ab8d36c9SPeter Brune PetscFunctionBegin; 762ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 763ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 764ab8d36c9SPeter Brune *mat = fas->interpolate; 765ab8d36c9SPeter Brune PetscFunctionReturn(0); 766ab8d36c9SPeter Brune } 767ab8d36c9SPeter Brune 768ab8d36c9SPeter Brune #undef __FUNCT__ 769ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 770ab8d36c9SPeter Brune /*@ 771ab8d36c9SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 772ab8d36c9SPeter Brune from level l to l-1. 773ab8d36c9SPeter Brune 774ab8d36c9SPeter Brune Input Parameters: 775ab8d36c9SPeter Brune + snes - the multigrid context 776ab8d36c9SPeter Brune . mat - the restriction matrix 777ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 778ab8d36c9SPeter Brune 779ab8d36c9SPeter Brune Level: advanced 780ab8d36c9SPeter Brune 781ab8d36c9SPeter Brune Notes: 782ab8d36c9SPeter Brune Usually this is the same matrix used also to set the interpolation 783ab8d36c9SPeter Brune for the same level. 784ab8d36c9SPeter Brune 785ab8d36c9SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 786ab8d36c9SPeter Brune out from the matrix size which one it is. 787ab8d36c9SPeter Brune 788ab8d36c9SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 789ab8d36c9SPeter Brune is used. 790ab8d36c9SPeter Brune 791ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 792ab8d36c9SPeter Brune 793ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 794ab8d36c9SPeter Brune @*/ 795ab8d36c9SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 796ab8d36c9SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 797ab8d36c9SPeter Brune PetscErrorCode ierr; 798ab8d36c9SPeter Brune SNES levelsnes; 799ab8d36c9SPeter Brune PetscFunctionBegin; 800ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 801ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 802ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 803ab8d36c9SPeter Brune ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 804ab8d36c9SPeter Brune fas->restrct = mat; 805ab8d36c9SPeter Brune PetscFunctionReturn(0); 806ab8d36c9SPeter Brune } 807ab8d36c9SPeter Brune 808ab8d36c9SPeter Brune #undef __FUNCT__ 809ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetRestriction" 810ab8d36c9SPeter Brune /*@ 811ab8d36c9SPeter Brune SNESFASGetRestriction - Gets the matrix used to calculate the 812ab8d36c9SPeter Brune restriction from l to the l-1th level 813ab8d36c9SPeter Brune 814ab8d36c9SPeter Brune Input Parameters: 815ab8d36c9SPeter Brune + snes - the multigrid context 816ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 817ab8d36c9SPeter Brune 818ab8d36c9SPeter Brune Output Parameters: 819ab8d36c9SPeter Brune . mat - the interpolation operator 820ab8d36c9SPeter Brune 821ab8d36c9SPeter Brune Level: advanced 822ab8d36c9SPeter Brune 823ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, restrict, level 824ab8d36c9SPeter Brune 825ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale() 826ab8d36c9SPeter Brune @*/ 827ab8d36c9SPeter Brune PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) { 828ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 829ab8d36c9SPeter Brune PetscErrorCode ierr; 830ab8d36c9SPeter Brune SNES levelsnes; 831ab8d36c9SPeter Brune PetscFunctionBegin; 832ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 833ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 834ab8d36c9SPeter Brune *mat = fas->restrct; 835ab8d36c9SPeter Brune PetscFunctionReturn(0); 836ab8d36c9SPeter Brune } 837ab8d36c9SPeter Brune 838ab8d36c9SPeter Brune 839ab8d36c9SPeter Brune #undef __FUNCT__ 840ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInjection" 841ab8d36c9SPeter Brune /*@ 842ab8d36c9SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 843ab8d36c9SPeter Brune from level l to l-1. 844ab8d36c9SPeter Brune 845ab8d36c9SPeter Brune Input Parameters: 846ab8d36c9SPeter Brune + snes - the multigrid context 847ab8d36c9SPeter Brune . mat - the restriction matrix 848ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 849ab8d36c9SPeter Brune 850ab8d36c9SPeter Brune Level: advanced 851ab8d36c9SPeter Brune 852ab8d36c9SPeter Brune Notes: 853ab8d36c9SPeter Brune If you do not set this, the restriction and rscale is used to 854ab8d36c9SPeter Brune project the solution instead. 855ab8d36c9SPeter Brune 856ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 857ab8d36c9SPeter Brune 858ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 859ab8d36c9SPeter Brune @*/ 860ab8d36c9SPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 861ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 862ab8d36c9SPeter Brune PetscErrorCode ierr; 863ab8d36c9SPeter Brune SNES levelsnes; 864ab8d36c9SPeter Brune PetscFunctionBegin; 865ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 866ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 867ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 868ab8d36c9SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 869ab8d36c9SPeter Brune fas->inject = mat; 870ab8d36c9SPeter Brune PetscFunctionReturn(0); 871ab8d36c9SPeter Brune } 872ab8d36c9SPeter Brune 873ab8d36c9SPeter Brune 874ab8d36c9SPeter Brune #undef __FUNCT__ 875ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInjection" 876ab8d36c9SPeter Brune /*@ 877ab8d36c9SPeter Brune SNESFASGetInjection - Gets the matrix used to calculate the 878ab8d36c9SPeter Brune injection from l-1 to the lth level 879ab8d36c9SPeter Brune 880ab8d36c9SPeter Brune Input Parameters: 881ab8d36c9SPeter Brune + snes - the multigrid context 882ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 883ab8d36c9SPeter Brune 884ab8d36c9SPeter Brune Output Parameters: 885ab8d36c9SPeter Brune . mat - the injection operator 886ab8d36c9SPeter Brune 887ab8d36c9SPeter Brune Level: advanced 888ab8d36c9SPeter Brune 889ab8d36c9SPeter Brune .keywords: FAS, multigrid, get, restrict, level 890ab8d36c9SPeter Brune 891ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale() 892ab8d36c9SPeter Brune @*/ 893ab8d36c9SPeter Brune PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) { 894ab8d36c9SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 895ab8d36c9SPeter Brune PetscErrorCode ierr; 896ab8d36c9SPeter Brune SNES levelsnes; 897ab8d36c9SPeter Brune PetscFunctionBegin; 898ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 899ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 900ab8d36c9SPeter Brune *mat = fas->inject; 901ab8d36c9SPeter Brune PetscFunctionReturn(0); 902ab8d36c9SPeter Brune } 903ab8d36c9SPeter Brune 904ab8d36c9SPeter Brune #undef __FUNCT__ 905ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 906ab8d36c9SPeter Brune /*@ 907ab8d36c9SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 908ab8d36c9SPeter Brune operator from level l to l-1. 909ab8d36c9SPeter Brune 910ab8d36c9SPeter Brune Input Parameters: 911ab8d36c9SPeter Brune + snes - the multigrid context 912ab8d36c9SPeter Brune . rscale - the restriction scaling 913ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 914ab8d36c9SPeter Brune 915ab8d36c9SPeter Brune Level: advanced 916ab8d36c9SPeter Brune 917ab8d36c9SPeter Brune Notes: 918ab8d36c9SPeter Brune This is only used in the case that the injection is not set. 919ab8d36c9SPeter Brune 920ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 921ab8d36c9SPeter Brune 922ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 923ab8d36c9SPeter Brune @*/ 924ab8d36c9SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 925ab8d36c9SPeter Brune SNES_FAS * fas; 926ab8d36c9SPeter Brune PetscErrorCode ierr; 927ab8d36c9SPeter Brune SNES levelsnes; 928ab8d36c9SPeter Brune PetscFunctionBegin; 929ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 930ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 931ab8d36c9SPeter Brune ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 932ab8d36c9SPeter Brune ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 933ab8d36c9SPeter Brune fas->rscale = rscale; 934ab8d36c9SPeter Brune PetscFunctionReturn(0); 935ab8d36c9SPeter Brune } 936ab8d36c9SPeter Brune 937ab8d36c9SPeter Brune #undef __FUNCT__ 938ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmoother" 939ab8d36c9SPeter Brune /*@ 940ab8d36c9SPeter Brune SNESFASGetSmoother - Gets the default smoother on a level. 941ab8d36c9SPeter Brune 942ab8d36c9SPeter Brune Input Parameters: 943ab8d36c9SPeter Brune + snes - the multigrid context 944ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 945ab8d36c9SPeter Brune 946ab8d36c9SPeter Brune Output Parameters: 947ab8d36c9SPeter Brune smooth - the smoother 948ab8d36c9SPeter Brune 949ab8d36c9SPeter Brune Level: advanced 950ab8d36c9SPeter Brune 951ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 952ab8d36c9SPeter Brune 953ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 954ab8d36c9SPeter Brune @*/ 955ab8d36c9SPeter Brune PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) { 956ab8d36c9SPeter Brune SNES_FAS * fas; 957ab8d36c9SPeter Brune PetscErrorCode ierr; 958ab8d36c9SPeter Brune SNES levelsnes; 959ab8d36c9SPeter Brune PetscFunctionBegin; 960ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 961ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 962ab8d36c9SPeter Brune if (!fas->smoothd) { 963ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 964ab8d36c9SPeter Brune } 965ab8d36c9SPeter Brune *smooth = fas->smoothd; 966ab8d36c9SPeter Brune PetscFunctionReturn(0); 967ab8d36c9SPeter Brune } 968ab8d36c9SPeter Brune 969ab8d36c9SPeter Brune #undef __FUNCT__ 970ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherDown" 971ab8d36c9SPeter Brune /*@ 972ab8d36c9SPeter Brune SNESFASGetSmootherDown - Gets the downsmoother on a level. 973ab8d36c9SPeter Brune 974ab8d36c9SPeter Brune Input Parameters: 975ab8d36c9SPeter Brune + snes - the multigrid context 976ab8d36c9SPeter Brune - level - the level (0 is coarsest) to supply 977ab8d36c9SPeter Brune 978ab8d36c9SPeter Brune Output Parameters: 979ab8d36c9SPeter Brune smooth - the smoother 980ab8d36c9SPeter Brune 981ab8d36c9SPeter Brune Level: advanced 982ab8d36c9SPeter Brune 983ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 984ab8d36c9SPeter Brune 985ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 986ab8d36c9SPeter Brune @*/ 987ab8d36c9SPeter Brune PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) { 988ab8d36c9SPeter Brune SNES_FAS * fas; 989ab8d36c9SPeter Brune PetscErrorCode ierr; 990ab8d36c9SPeter Brune SNES levelsnes; 991ab8d36c9SPeter Brune PetscFunctionBegin; 992ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 993ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 994ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 995ab8d36c9SPeter Brune if (!fas->smoothd) { 996ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 997ab8d36c9SPeter Brune } 998ab8d36c9SPeter Brune if (!fas->smoothu) { 999ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 1000ab8d36c9SPeter Brune } 1001ab8d36c9SPeter Brune *smooth = fas->smoothd; 1002ab8d36c9SPeter Brune PetscFunctionReturn(0); 1003ab8d36c9SPeter Brune } 1004ab8d36c9SPeter Brune 1005ab8d36c9SPeter Brune #undef __FUNCT__ 1006ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherUp" 1007ab8d36c9SPeter Brune /*@ 1008ab8d36c9SPeter Brune SNESFASGetSmootherUp - Gets the upsmoother on a level. 1009ab8d36c9SPeter Brune 1010ab8d36c9SPeter Brune Input Parameters: 1011ab8d36c9SPeter Brune + snes - the multigrid context 1012ab8d36c9SPeter Brune - level - the level (0 is coarsest) 1013ab8d36c9SPeter Brune 1014ab8d36c9SPeter Brune Output Parameters: 1015ab8d36c9SPeter Brune smooth - the smoother 1016ab8d36c9SPeter Brune 1017ab8d36c9SPeter Brune Level: advanced 1018ab8d36c9SPeter Brune 1019ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level 1020ab8d36c9SPeter Brune 1021ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1022ab8d36c9SPeter Brune @*/ 1023ab8d36c9SPeter Brune PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) { 1024ab8d36c9SPeter Brune SNES_FAS * fas; 1025ab8d36c9SPeter Brune PetscErrorCode ierr; 1026ab8d36c9SPeter Brune SNES levelsnes; 1027ab8d36c9SPeter Brune PetscFunctionBegin; 1028ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr); 1029ab8d36c9SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1030ab8d36c9SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1031ab8d36c9SPeter Brune if (!fas->smoothd) { 1032ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 1033ab8d36c9SPeter Brune } 1034ab8d36c9SPeter Brune if (!fas->smoothu) { 1035ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 1036ab8d36c9SPeter Brune } 1037ab8d36c9SPeter Brune *smooth = fas->smoothu; 1038ab8d36c9SPeter Brune PetscFunctionReturn(0); 1039ab8d36c9SPeter Brune } 1040d6ad1212SPeter Brune 1041d6ad1212SPeter Brune #undef __FUNCT__ 1042d6ad1212SPeter Brune #define __FUNCT__ "SNESFASGetCoarseSolve" 1043d6ad1212SPeter Brune /*@ 1044d6ad1212SPeter Brune SNESFASGetCoarseSolve - Gets the coarsest solver. 1045d6ad1212SPeter Brune 1046d6ad1212SPeter Brune Input Parameters: 1047d6ad1212SPeter Brune + snes - the multigrid context 1048d6ad1212SPeter Brune 1049d6ad1212SPeter Brune Output Parameters: 1050d6ad1212SPeter Brune solve - the coarse-level solver 1051d6ad1212SPeter Brune 1052d6ad1212SPeter Brune Level: advanced 1053d6ad1212SPeter Brune 1054d6ad1212SPeter Brune .keywords: FAS, MG, get, multigrid, solver, coarse 1055d6ad1212SPeter Brune 1056d6ad1212SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 1057d6ad1212SPeter Brune @*/ 1058d6ad1212SPeter Brune PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *smooth) { 1059d6ad1212SPeter Brune SNES_FAS * fas; 1060d6ad1212SPeter Brune PetscErrorCode ierr; 1061d6ad1212SPeter Brune SNES levelsnes; 1062d6ad1212SPeter Brune PetscFunctionBegin; 1063d6ad1212SPeter Brune ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr); 1064d6ad1212SPeter Brune fas = (SNES_FAS *)levelsnes->data; 1065d6ad1212SPeter Brune /* if the user chooses to differentiate smoothers, create them both at this point */ 1066d6ad1212SPeter Brune if (!fas->smoothd) { 1067d6ad1212SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 1068d6ad1212SPeter Brune } 1069d6ad1212SPeter Brune *smooth = fas->smoothd; 1070d6ad1212SPeter Brune PetscFunctionReturn(0); 1071d6ad1212SPeter Brune } 1072