1421d9b32SPeter Brune /* Defines the basic SNES object */ 26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3421d9b32SPeter Brune 4421d9b32SPeter Brune /*MC 5421d9b32SPeter Brune Full Approximation Scheme nonlinear multigrid solver. 6421d9b32SPeter Brune 7421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 8421d9b32SPeter Brune 9421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 10421d9b32SPeter Brune M*/ 11421d9b32SPeter Brune 12421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes); 13421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 14421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 15421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 16421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 17421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 186273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 19421d9b32SPeter Brune 20421d9b32SPeter Brune EXTERN_C_BEGIN 21421d9b32SPeter Brune 22421d9b32SPeter Brune #undef __FUNCT__ 23421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 24421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes) 25421d9b32SPeter Brune { 26421d9b32SPeter Brune SNES_FAS * fas; 27421d9b32SPeter Brune PetscErrorCode ierr; 28421d9b32SPeter Brune 29421d9b32SPeter Brune PetscFunctionBegin; 30421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 31421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 32421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 33421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 34421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 35421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 36421d9b32SPeter Brune 37ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 38ed020824SBarry Smith snes->usespc = PETSC_FALSE; 39ed020824SBarry Smith 40421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 41421d9b32SPeter Brune snes->data = (void*) fas; 42421d9b32SPeter Brune fas->level = 0; 43293a7e31SPeter Brune fas->levels = 1; 44ee78dd50SPeter Brune fas->n_cycles = 1; 45ee78dd50SPeter Brune fas->max_up_it = 1; 46ee78dd50SPeter Brune fas->max_down_it = 1; 47ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 48ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 49421d9b32SPeter Brune fas->next = PETSC_NULL; 506273346dSPeter Brune fas->previous = PETSC_NULL; 51421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 52421d9b32SPeter Brune fas->restrct = PETSC_NULL; 53efe1f98aSPeter Brune fas->inject = PETSC_NULL; 54cc05f883SPeter Brune fas->monitor = PETSC_NULL; 55cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 56421d9b32SPeter Brune PetscFunctionReturn(0); 57421d9b32SPeter Brune } 58421d9b32SPeter Brune EXTERN_C_END 59421d9b32SPeter Brune 60421d9b32SPeter Brune #undef __FUNCT__ 61421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 62c79ef259SPeter Brune /*@ 63722262beSPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS. 64c79ef259SPeter Brune 65c79ef259SPeter Brune Input Parameter: 66c79ef259SPeter Brune . snes - the preconditioner context 67c79ef259SPeter Brune 68c79ef259SPeter Brune Output parameter: 69c79ef259SPeter Brune . levels - the number of levels 70c79ef259SPeter Brune 71c79ef259SPeter Brune Level: advanced 72c79ef259SPeter Brune 73c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 74c79ef259SPeter Brune 75c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 76c79ef259SPeter Brune @*/ 77421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 78421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 79421d9b32SPeter Brune PetscFunctionBegin; 80ee1fd11aSPeter Brune *levels = fas->levels; 81421d9b32SPeter Brune PetscFunctionReturn(0); 82421d9b32SPeter Brune } 83421d9b32SPeter Brune 84421d9b32SPeter Brune #undef __FUNCT__ 85646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 86c79ef259SPeter Brune /*@ 87c79ef259SPeter Brune SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 88c79ef259SPeter Brune complicated cycling. 89c79ef259SPeter Brune 90c79ef259SPeter Brune Logically Collective on SNES 91c79ef259SPeter Brune 92c79ef259SPeter Brune Input Parameters: 93c79ef259SPeter Brune + snes - the multigrid context 94c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 95c79ef259SPeter Brune 96c79ef259SPeter Brune Options Database Key: 97c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 98c79ef259SPeter Brune 99c79ef259SPeter Brune Level: advanced 100c79ef259SPeter Brune 101c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 102c79ef259SPeter Brune 103c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 104c79ef259SPeter Brune @*/ 105646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 106646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 107646217ecSPeter Brune PetscErrorCode ierr; 108646217ecSPeter Brune PetscFunctionBegin; 109646217ecSPeter Brune fas->n_cycles = cycles; 110646217ecSPeter Brune if (fas->next) { 111646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 112646217ecSPeter Brune } 113646217ecSPeter Brune PetscFunctionReturn(0); 114646217ecSPeter Brune } 115646217ecSPeter Brune 116eff52c0eSPeter Brune #undef __FUNCT__ 117c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 118c79ef259SPeter Brune /*@ 119722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 120c79ef259SPeter Brune 121c79ef259SPeter Brune Logically Collective on SNES 122c79ef259SPeter Brune 123c79ef259SPeter Brune Input Parameters: 124c79ef259SPeter Brune + snes - the multigrid context 125c79ef259SPeter Brune . level - the level to set the number of cycles on 126c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 127c79ef259SPeter Brune 128c79ef259SPeter Brune Level: advanced 129c79ef259SPeter Brune 130c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 131c79ef259SPeter Brune 132c79ef259SPeter Brune .seealso: SNESFASSetCycles() 133c79ef259SPeter Brune @*/ 134c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 135c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 136c79ef259SPeter Brune PetscInt top_level = fas->level,i; 137c79ef259SPeter Brune 138c79ef259SPeter Brune PetscFunctionBegin; 139c79ef259SPeter Brune if (level > top_level) 140c79ef259SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 141c79ef259SPeter Brune /* get to the correct level */ 142c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 143c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 144c79ef259SPeter Brune } 145c79ef259SPeter Brune if (fas->level != level) 146c79ef259SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 147c79ef259SPeter Brune fas->n_cycles = cycles; 148c79ef259SPeter Brune PetscFunctionReturn(0); 149c79ef259SPeter Brune } 150c79ef259SPeter Brune 151c79ef259SPeter Brune #undef __FUNCT__ 152eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 153aeed3662SMatthew G Knepley /*@C 154c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 155c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 156c79ef259SPeter Brune and nonlinear preconditioners. 157c79ef259SPeter Brune 158c79ef259SPeter Brune Logically Collective on SNES 159c79ef259SPeter Brune 160c79ef259SPeter Brune Input Parameters: 161c79ef259SPeter Brune + snes - the multigrid context 162c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 163c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 164c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 165c79ef259SPeter Brune 166c79ef259SPeter Brune Level: advanced 167c79ef259SPeter Brune 168c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 169c79ef259SPeter Brune 170c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 171c79ef259SPeter Brune @*/ 172eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 173eff52c0eSPeter Brune PetscErrorCode ierr = 0; 174eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 175eff52c0eSPeter Brune PetscFunctionBegin; 176eff52c0eSPeter Brune 177eff52c0eSPeter Brune /* use or don't use it according to user wishes*/ 178d28543b3SPeter Brune snes->usegs = use_gs; 179eff52c0eSPeter Brune if (gsfunc) { 180eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 181eff52c0eSPeter Brune /* push the provided GS up the tree */ 182eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 183eff52c0eSPeter Brune } else if (snes->ops->computegs) { 184eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 185eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 186eff52c0eSPeter Brune } else if (use_gs) { 187eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 188eff52c0eSPeter Brune } 189eff52c0eSPeter Brune PetscFunctionReturn(0); 190eff52c0eSPeter Brune } 191eff52c0eSPeter Brune 192eff52c0eSPeter Brune #undef __FUNCT__ 193eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 194aeed3662SMatthew G Knepley /*@C 195c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 196c79ef259SPeter Brune 197c79ef259SPeter Brune Logically Collective on SNES 198c79ef259SPeter Brune 199c79ef259SPeter Brune Input Parameters: 200c79ef259SPeter Brune + snes - the multigrid context 201c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 202c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 203c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 204c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 205c79ef259SPeter Brune 206c79ef259SPeter Brune Level: advanced 207c79ef259SPeter Brune 208c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 209c79ef259SPeter Brune 210c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 211c79ef259SPeter Brune @*/ 212eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 213eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 214eff52c0eSPeter Brune PetscErrorCode ierr; 215eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 216eff52c0eSPeter Brune SNES cur_snes = snes; 217eff52c0eSPeter Brune PetscFunctionBegin; 218eff52c0eSPeter Brune if (level > top_level) 219eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 220eff52c0eSPeter Brune /* get to the correct level */ 221eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 222eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 223eff52c0eSPeter Brune cur_snes = fas->next; 224eff52c0eSPeter Brune } 225eff52c0eSPeter Brune if (fas->level != level) 226eff52c0eSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 227d28543b3SPeter Brune snes->usegs = use_gs; 228eff52c0eSPeter Brune if (gsfunc) { 2296273346dSPeter Brune ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 230eff52c0eSPeter Brune } 231eff52c0eSPeter Brune PetscFunctionReturn(0); 232eff52c0eSPeter Brune } 233eff52c0eSPeter Brune 234646217ecSPeter Brune #undef __FUNCT__ 235421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 236c79ef259SPeter Brune /*@ 237c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 238c79ef259SPeter Brune level of the FAS hierarchy. 239c79ef259SPeter Brune 240c79ef259SPeter Brune Input Parameters: 241c79ef259SPeter Brune + snes - the multigrid context 242c79ef259SPeter Brune level - the level to get 243c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 244c79ef259SPeter Brune 245c79ef259SPeter Brune Level: advanced 246c79ef259SPeter Brune 247c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 248c79ef259SPeter Brune 249c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 250c79ef259SPeter Brune @*/ 251421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 252421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 253421d9b32SPeter Brune PetscInt levels = fas->level; 254421d9b32SPeter Brune PetscInt i; 255421d9b32SPeter Brune PetscFunctionBegin; 256421d9b32SPeter Brune *lsnes = snes; 257421d9b32SPeter Brune if (fas->level < level) { 258421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 259421d9b32SPeter Brune } 260421d9b32SPeter Brune if (level > levels - 1) { 261421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 262421d9b32SPeter Brune } 263421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 264421d9b32SPeter Brune *lsnes = fas->next; 265421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 266421d9b32SPeter Brune } 267421d9b32SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 268421d9b32SPeter Brune PetscFunctionReturn(0); 269421d9b32SPeter Brune } 270421d9b32SPeter Brune 271421d9b32SPeter Brune #undef __FUNCT__ 272421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 273aeed3662SMatthew G Knepley /*@C 274c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 275c79ef259SPeter Brune Must be called before any other FAS routine. 276c79ef259SPeter Brune 277c79ef259SPeter Brune Input Parameters: 278c79ef259SPeter Brune + snes - the snes context 279c79ef259SPeter Brune . levels - the number of levels 280c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 281c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 282c79ef259SPeter Brune Fortran. 283c79ef259SPeter Brune 284c79ef259SPeter Brune Level: intermediate 285c79ef259SPeter Brune 286c79ef259SPeter Brune Notes: 287c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 288c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 289c79ef259SPeter Brune 290c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 291c79ef259SPeter Brune 292c79ef259SPeter Brune .seealso: SNESFASGetLevels() 293c79ef259SPeter Brune @*/ 294421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 295421d9b32SPeter Brune PetscErrorCode ierr; 296421d9b32SPeter Brune PetscInt i; 297421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 2986273346dSPeter Brune SNES prevsnes; 299421d9b32SPeter Brune MPI_Comm comm; 300421d9b32SPeter Brune PetscFunctionBegin; 301ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 302ee1fd11aSPeter Brune if (levels == fas->levels) { 303ee1fd11aSPeter Brune if (!comms) 304ee1fd11aSPeter Brune PetscFunctionReturn(0); 305ee1fd11aSPeter Brune } 306421d9b32SPeter Brune /* user has changed the number of levels; reset */ 307421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 308421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 309421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 310ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3116273346dSPeter Brune fas->previous = PETSC_NULL; 3126273346dSPeter Brune prevsnes = snes; 313421d9b32SPeter Brune /* setup the finest level */ 314421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 315421d9b32SPeter Brune if (comms) comm = comms[i]; 316421d9b32SPeter Brune fas->level = i; 317421d9b32SPeter Brune fas->levels = levels; 318ee1fd11aSPeter Brune fas->next = PETSC_NULL; 319421d9b32SPeter Brune if (i > 0) { 320421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3214a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 322421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 323794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3246273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3256273346dSPeter Brune prevsnes = fas->next; 3266273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 327421d9b32SPeter Brune } 328421d9b32SPeter Brune } 329421d9b32SPeter Brune PetscFunctionReturn(0); 330421d9b32SPeter Brune } 331421d9b32SPeter Brune 332421d9b32SPeter Brune #undef __FUNCT__ 333c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 334c79ef259SPeter Brune /*@ 335c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 336c79ef259SPeter Brune use on all levels. 337c79ef259SPeter Brune 338c79ef259SPeter Brune Logically Collective on PC 339c79ef259SPeter Brune 340c79ef259SPeter Brune Input Parameters: 341c79ef259SPeter Brune + snes - the multigrid context 342c79ef259SPeter Brune - n - the number of smoothing steps 343c79ef259SPeter Brune 344c79ef259SPeter Brune Options Database Key: 345d28543b3SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 346c79ef259SPeter Brune 347c79ef259SPeter Brune Level: advanced 348c79ef259SPeter Brune 349c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 350c79ef259SPeter Brune 351c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 352c79ef259SPeter Brune @*/ 353c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 354c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 355c79ef259SPeter Brune PetscErrorCode ierr = 0; 356c79ef259SPeter Brune PetscFunctionBegin; 357c79ef259SPeter Brune fas->max_up_it = n; 358c79ef259SPeter Brune if (fas->next) { 359c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 360c79ef259SPeter Brune } 361c79ef259SPeter Brune PetscFunctionReturn(0); 362c79ef259SPeter Brune } 363c79ef259SPeter Brune 364c79ef259SPeter Brune #undef __FUNCT__ 365c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 366c79ef259SPeter Brune /*@ 367c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 368c79ef259SPeter Brune use on all levels. 369c79ef259SPeter Brune 370c79ef259SPeter Brune Logically Collective on PC 371c79ef259SPeter Brune 372c79ef259SPeter Brune Input Parameters: 373c79ef259SPeter Brune + snes - the multigrid context 374c79ef259SPeter Brune - n - the number of smoothing steps 375c79ef259SPeter Brune 376c79ef259SPeter Brune Options Database Key: 377d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 378c79ef259SPeter Brune 379c79ef259SPeter Brune Level: advanced 380c79ef259SPeter Brune 381c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 382c79ef259SPeter Brune 383c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 384c79ef259SPeter Brune @*/ 385c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 386c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 387c79ef259SPeter Brune PetscErrorCode ierr = 0; 388c79ef259SPeter Brune PetscFunctionBegin; 389c79ef259SPeter Brune fas->max_down_it = n; 390c79ef259SPeter Brune if (fas->next) { 391c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 392c79ef259SPeter Brune } 393c79ef259SPeter Brune PetscFunctionReturn(0); 394c79ef259SPeter Brune } 395c79ef259SPeter Brune 396c79ef259SPeter Brune #undef __FUNCT__ 397421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 398c79ef259SPeter Brune /*@ 399c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 400c79ef259SPeter Brune interpolation from l-1 to the lth level 401c79ef259SPeter Brune 402c79ef259SPeter Brune Input Parameters: 403c79ef259SPeter Brune + snes - the multigrid context 404c79ef259SPeter Brune . mat - the interpolation operator 405c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 406c79ef259SPeter Brune 407c79ef259SPeter Brune Level: advanced 408c79ef259SPeter Brune 409c79ef259SPeter Brune Notes: 410c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 411c79ef259SPeter Brune for the same level. 412c79ef259SPeter Brune 413c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 414c79ef259SPeter Brune out from the matrix size which one it is. 415c79ef259SPeter Brune 416c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 417c79ef259SPeter Brune 418c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale() 419c79ef259SPeter Brune @*/ 420421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 421421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 422d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 423421d9b32SPeter Brune 424421d9b32SPeter Brune PetscFunctionBegin; 425421d9b32SPeter Brune if (level > top_level) 426421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 427421d9b32SPeter Brune /* get to the correct level */ 428d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 429421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 430421d9b32SPeter Brune } 431421d9b32SPeter Brune if (fas->level != level) 432421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 433421d9b32SPeter Brune fas->interpolate = mat; 434421d9b32SPeter Brune PetscFunctionReturn(0); 435421d9b32SPeter Brune } 436421d9b32SPeter Brune 437421d9b32SPeter Brune #undef __FUNCT__ 438421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 439c79ef259SPeter Brune /*@ 440c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 441c79ef259SPeter Brune from level l to l-1. 442c79ef259SPeter Brune 443c79ef259SPeter Brune Input Parameters: 444c79ef259SPeter Brune + snes - the multigrid context 445c79ef259SPeter Brune . mat - the restriction matrix 446c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 447c79ef259SPeter Brune 448c79ef259SPeter Brune Level: advanced 449c79ef259SPeter Brune 450c79ef259SPeter Brune Notes: 451c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 452c79ef259SPeter Brune for the same level. 453c79ef259SPeter Brune 454c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 455c79ef259SPeter Brune out from the matrix size which one it is. 456c79ef259SPeter Brune 457c79ef259SPeter Brune If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 458c79ef259SPeter Brune is used. 459c79ef259SPeter Brune 460c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 461c79ef259SPeter Brune 462c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 463c79ef259SPeter Brune @*/ 464421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 465421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 466d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 467421d9b32SPeter Brune 468421d9b32SPeter Brune PetscFunctionBegin; 469421d9b32SPeter Brune if (level > top_level) 470421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 471421d9b32SPeter Brune /* get to the correct level */ 472d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 473421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 474421d9b32SPeter Brune } 475421d9b32SPeter Brune if (fas->level != level) 476421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 477421d9b32SPeter Brune fas->restrct = mat; 478421d9b32SPeter Brune PetscFunctionReturn(0); 479421d9b32SPeter Brune } 480421d9b32SPeter Brune 481efe1f98aSPeter Brune 482421d9b32SPeter Brune #undef __FUNCT__ 483421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 484c79ef259SPeter Brune /*@ 485c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 486c79ef259SPeter Brune operator from level l to l-1. 487c79ef259SPeter Brune 488c79ef259SPeter Brune Input Parameters: 489c79ef259SPeter Brune + snes - the multigrid context 490c79ef259SPeter Brune . rscale - the restriction scaling 491c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 492c79ef259SPeter Brune 493c79ef259SPeter Brune Level: advanced 494c79ef259SPeter Brune 495c79ef259SPeter Brune Notes: 496c79ef259SPeter Brune This is only used in the case that the injection is not set. 497c79ef259SPeter Brune 498c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 499c79ef259SPeter Brune 500c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 501c79ef259SPeter Brune @*/ 502421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 503421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 504d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 505421d9b32SPeter Brune 506421d9b32SPeter Brune PetscFunctionBegin; 507421d9b32SPeter Brune if (level > top_level) 508421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 509421d9b32SPeter Brune /* get to the correct level */ 510d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 511421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 512421d9b32SPeter Brune } 513421d9b32SPeter Brune if (fas->level != level) 514421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 515421d9b32SPeter Brune fas->rscale = rscale; 516421d9b32SPeter Brune PetscFunctionReturn(0); 517421d9b32SPeter Brune } 518421d9b32SPeter Brune 519efe1f98aSPeter Brune 520efe1f98aSPeter Brune #undef __FUNCT__ 521efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 522c79ef259SPeter Brune /*@ 523c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 524c79ef259SPeter Brune from level l to l-1. 525c79ef259SPeter Brune 526c79ef259SPeter Brune Input Parameters: 527c79ef259SPeter Brune + snes - the multigrid context 528c79ef259SPeter Brune . mat - the restriction matrix 529c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 530c79ef259SPeter Brune 531c79ef259SPeter Brune Level: advanced 532c79ef259SPeter Brune 533c79ef259SPeter Brune Notes: 534c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 535c79ef259SPeter Brune project the solution instead. 536c79ef259SPeter Brune 537c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 538c79ef259SPeter Brune 539c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 540c79ef259SPeter Brune @*/ 541efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 542efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 543efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 544efe1f98aSPeter Brune 545efe1f98aSPeter Brune PetscFunctionBegin; 546efe1f98aSPeter Brune if (level > top_level) 547efe1f98aSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 548efe1f98aSPeter Brune /* get to the correct level */ 549efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 550efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 551efe1f98aSPeter Brune } 552efe1f98aSPeter Brune if (fas->level != level) 553efe1f98aSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 554efe1f98aSPeter Brune fas->inject = mat; 555efe1f98aSPeter Brune PetscFunctionReturn(0); 556efe1f98aSPeter Brune } 557efe1f98aSPeter Brune 558421d9b32SPeter Brune #undef __FUNCT__ 559421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 560421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 561421d9b32SPeter Brune { 56277df8cc4SPeter Brune PetscErrorCode ierr = 0; 563421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 564421d9b32SPeter Brune 565421d9b32SPeter Brune PetscFunctionBegin; 566742fe5e2SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 567742fe5e2SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 568742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 569421d9b32SPeter Brune PetscFunctionReturn(0); 570421d9b32SPeter Brune } 571421d9b32SPeter Brune 572421d9b32SPeter Brune #undef __FUNCT__ 573421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 574421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 575421d9b32SPeter Brune { 576421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 577742fe5e2SPeter Brune PetscErrorCode ierr = 0; 578421d9b32SPeter Brune 579421d9b32SPeter Brune PetscFunctionBegin; 580421d9b32SPeter Brune /* recursively resets and then destroys */ 58151e86f29SPeter Brune if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 58251e86f29SPeter Brune if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 583efe1f98aSPeter Brune if (fas->inject) { 584efe1f98aSPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 585efe1f98aSPeter Brune } 58651e86f29SPeter Brune if (fas->interpolate == fas->restrct) { 58751e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 58851e86f29SPeter Brune fas->restrct = PETSC_NULL; 58951e86f29SPeter Brune } else { 59051e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 59151e86f29SPeter Brune if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 59251e86f29SPeter Brune } 59351e86f29SPeter Brune if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 594421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 595421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 596ee1fd11aSPeter Brune 597421d9b32SPeter Brune PetscFunctionReturn(0); 598421d9b32SPeter Brune } 599421d9b32SPeter Brune 600421d9b32SPeter Brune #undef __FUNCT__ 601421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 602421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 603421d9b32SPeter Brune { 60448bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 605421d9b32SPeter Brune PetscErrorCode ierr; 606efe1f98aSPeter Brune VecScatter injscatter; 607d1adcc6fSPeter Brune PetscInt dm_levels; 608d1adcc6fSPeter Brune 609421d9b32SPeter Brune PetscFunctionBegin; 610eff52c0eSPeter Brune 611cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 612d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 613d1adcc6fSPeter Brune dm_levels++; 614cc05f883SPeter Brune if (dm_levels > fas->levels) { 615d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 616cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 617d1adcc6fSPeter Brune } 618d1adcc6fSPeter Brune } 619d1adcc6fSPeter Brune 6202f7ea302SPeter Brune ierr = SNESDefaultGetWork(snes, 1);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 621cc05f883SPeter Brune 62248bfdf8aSPeter Brune /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 62348bfdf8aSPeter Brune if (fas->upsmooth) { 62448bfdf8aSPeter Brune ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 62548bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 62648bfdf8aSPeter Brune ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 62748bfdf8aSPeter Brune } 62848bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 62948bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 63048bfdf8aSPeter Brune } 63148bfdf8aSPeter Brune if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 63248bfdf8aSPeter Brune ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 63348bfdf8aSPeter Brune } 63448bfdf8aSPeter Brune } 63548bfdf8aSPeter Brune if (fas->downsmooth) { 63648bfdf8aSPeter Brune ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 63748bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 63848bfdf8aSPeter Brune ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 63948bfdf8aSPeter Brune } 64048bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 64148bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 64248bfdf8aSPeter Brune } 64348bfdf8aSPeter Brune if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 64448bfdf8aSPeter Brune ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 64548bfdf8aSPeter Brune } 6461a266240SBarry Smith } 64748bfdf8aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 648646217ecSPeter Brune if (fas->next) { 6496273346dSPeter Brune if (fas->galerkin) { 6506273346dSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 6516273346dSPeter Brune } else { 65248bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->next->ops->computefunction) { 65348bfdf8aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 65448bfdf8aSPeter Brune } 65548bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 65648bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 65748bfdf8aSPeter Brune } 65848bfdf8aSPeter Brune if (snes->ops->computegs && !fas->next->ops->computegs) { 659646217ecSPeter Brune ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 660646217ecSPeter Brune } 661646217ecSPeter Brune } 6626273346dSPeter Brune } 663421d9b32SPeter Brune if (snes->dm) { 664421d9b32SPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 665421d9b32SPeter Brune if (fas->next) { 666421d9b32SPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 667ee78dd50SPeter Brune if (!fas->next->dm) { 668ee78dd50SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 669ee78dd50SPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 670ee78dd50SPeter Brune } 671421d9b32SPeter Brune /* set the interpolation and restriction from the DM */ 672ee78dd50SPeter Brune if (!fas->interpolate) { 673e727c939SJed Brown ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 674421d9b32SPeter Brune fas->restrct = fas->interpolate; 675421d9b32SPeter Brune } 676efe1f98aSPeter Brune /* set the injection from the DM */ 677efe1f98aSPeter Brune if (!fas->inject) { 678e727c939SJed Brown ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 679efe1f98aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 680efe1f98aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 681efe1f98aSPeter Brune } 682421d9b32SPeter Brune } 683ee78dd50SPeter Brune /* set the DMs of the pre and post-smoothers here */ 684*6bed07a3SBarry Smith if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 685*6bed07a3SBarry Smith if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 686421d9b32SPeter Brune } 687cc05f883SPeter Brune 6886273346dSPeter Brune /* setup FAS work vectors */ 6896273346dSPeter Brune if (fas->galerkin) { 6906273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 6916273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 6926273346dSPeter Brune } 6936273346dSPeter Brune 694fa9694d7SPeter Brune if (fas->next) { 695fa9694d7SPeter Brune /* gotta set up the solution vector for this to work */ 696*6bed07a3SBarry Smith if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 697742fe5e2SPeter Brune if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);} 698fa9694d7SPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 699fa9694d7SPeter Brune } 700fa9694d7SPeter Brune /* got to set them all up at once */ 701421d9b32SPeter Brune PetscFunctionReturn(0); 702421d9b32SPeter Brune } 703421d9b32SPeter Brune 704421d9b32SPeter Brune #undef __FUNCT__ 705421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 706421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 707421d9b32SPeter Brune { 708ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 709ee78dd50SPeter Brune PetscInt levels = 1; 710d1adcc6fSPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg; 711421d9b32SPeter Brune PetscErrorCode ierr; 712ee78dd50SPeter Brune const char *def_smooth = SNESNRICHARDSON; 713ee78dd50SPeter Brune char pre_type[256]; 714ee78dd50SPeter Brune char post_type[256]; 715ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 716421d9b32SPeter Brune 717421d9b32SPeter Brune PetscFunctionBegin; 718c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 719ee78dd50SPeter Brune 720ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 721ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 722ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 723c732cbdbSBarry Smith if (!flg && snes->dm) { 724c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 725c732cbdbSBarry Smith levels++; 726d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 727c732cbdbSBarry Smith } 728ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 729ee78dd50SPeter Brune } 730ee78dd50SPeter Brune 731ee78dd50SPeter Brune /* type of pre and/or post smoothers -- set both at once */ 732ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 733ee78dd50SPeter Brune ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 734d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 735d1adcc6fSPeter Brune if (smoothflg) { 736ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 737ee78dd50SPeter Brune } else { 738d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr); 739d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr); 740ee78dd50SPeter Brune } 741ee78dd50SPeter Brune 742ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 743d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 744d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 745d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 746ee78dd50SPeter Brune 7476273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 7486273346dSPeter Brune 749646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 750ee78dd50SPeter Brune 751ee78dd50SPeter Brune /* other options for the coarsest level */ 752ee78dd50SPeter Brune if (fas->level == 0) { 753d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 754ee78dd50SPeter Brune } 755ee78dd50SPeter Brune 756421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 7578cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 758eff52c0eSPeter Brune 759d28543b3SPeter Brune if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->usegs)) { 7608cc86e31SPeter Brune const char *prefix; 7618cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 7628cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 7638cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 7648cc86e31SPeter Brune if (fas->level || (fas->levels == 1)) { 765eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr); 7668cc86e31SPeter Brune } else { 7678cc86e31SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 7688cc86e31SPeter Brune } 7698cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 7708cc86e31SPeter Brune ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 7718cc86e31SPeter Brune } 7728cc86e31SPeter Brune 773d28543b3SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->usegs)) { 77467339d5cSBarry Smith const char *prefix; 77567339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 776ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 77767339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 778eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr); 779293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 780ee78dd50SPeter Brune ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 781ee78dd50SPeter Brune } 7821a266240SBarry Smith if (fas->upsmooth) { 7831a266240SBarry Smith ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 7841a266240SBarry Smith } 7851a266240SBarry Smith 7861a266240SBarry Smith if (fas->downsmooth) { 7871a266240SBarry Smith ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 7881a266240SBarry Smith } 789ee78dd50SPeter Brune 790742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 791742fe5e2SPeter Brune ierr = SNESSetTolerances(snes, 0.0, 0.0, 0.0, fas->n_cycles, 1000);CHKERRQ(ierr); 792742fe5e2SPeter Brune } 793742fe5e2SPeter Brune 794ee78dd50SPeter Brune if (monflg) { 795646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 796794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 7972f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 798742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 799293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 800293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 801d28543b3SPeter Brune } else { 802d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 803d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 804d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 805d28543b3SPeter Brune } 806ee78dd50SPeter Brune } 807ee78dd50SPeter Brune 808ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 809d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 810421d9b32SPeter Brune PetscFunctionReturn(0); 811421d9b32SPeter Brune } 812421d9b32SPeter Brune 813421d9b32SPeter Brune #undef __FUNCT__ 814421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 815421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 816421d9b32SPeter Brune { 817421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 818421d9b32SPeter Brune PetscBool iascii; 819421d9b32SPeter Brune PetscErrorCode ierr; 820421d9b32SPeter Brune PetscInt levels = fas->levels; 821421d9b32SPeter Brune PetscInt i; 822421d9b32SPeter Brune 823421d9b32SPeter Brune PetscFunctionBegin; 824421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 825421d9b32SPeter Brune if (iascii) { 826421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 827421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 828421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 829421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 830ee78dd50SPeter Brune if (fas->upsmooth) { 831ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 832421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 833ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 834421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 835421d9b32SPeter Brune } else { 836ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 837421d9b32SPeter Brune } 838ee78dd50SPeter Brune if (fas->downsmooth) { 839ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 840421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 841ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 842421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 843421d9b32SPeter Brune } else { 844ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 845421d9b32SPeter Brune } 846d28543b3SPeter Brune if (snes->usegs) { 847eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n", fas->level);CHKERRQ(ierr); 848eff52c0eSPeter Brune } 849421d9b32SPeter Brune if (fas->next) fas = (SNES_FAS *)fas->next->data; 850421d9b32SPeter Brune } 851421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 852421d9b32SPeter Brune } else { 853421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 854421d9b32SPeter Brune } 855421d9b32SPeter Brune PetscFunctionReturn(0); 856421d9b32SPeter Brune } 857421d9b32SPeter Brune 858421d9b32SPeter Brune /* 859fa9694d7SPeter Brune 860fa9694d7SPeter Brune Defines the FAS cycle as: 861fa9694d7SPeter Brune 862fa9694d7SPeter Brune fine problem: F(x) = 0 863fa9694d7SPeter Brune coarse problem: F^c(x) = b^c 864fa9694d7SPeter Brune 865fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 866fa9694d7SPeter Brune 867fa9694d7SPeter Brune correction: 868fa9694d7SPeter Brune 869fa9694d7SPeter Brune x = x + I(x^c - Rx) 870fa9694d7SPeter Brune 871421d9b32SPeter Brune */ 872fa9694d7SPeter Brune 873421d9b32SPeter Brune #undef __FUNCT__ 874421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private" 875646217ecSPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec X) { 876fa9694d7SPeter Brune 877fa9694d7SPeter Brune PetscErrorCode ierr; 878646217ecSPeter Brune Vec X_c, Xo_c, F_c, B_c,F,B; 879fa9694d7SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 880742fe5e2SPeter Brune PetscInt k; 8812d15c84fSPeter Brune PetscReal fnorm; 882742fe5e2SPeter Brune SNESConvergedReason reason; 883421d9b32SPeter Brune 884421d9b32SPeter Brune PetscFunctionBegin; 885f5a6d4f9SBarry Smith F = snes->vec_func; 886646217ecSPeter Brune B = snes->vec_rhs; 887fa9694d7SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 888d1adcc6fSPeter Brune if (fas->downsmooth) { 889d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 890742fe5e2SPeter Brune /* check convergence reason for the smoother */ 891742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 892742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 893742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 894742fe5e2SPeter Brune PetscFunctionReturn(0); 895742fe5e2SPeter Brune } 896d28543b3SPeter Brune } else if (snes->usegs && snes->ops->computegs) { 897794bee33SPeter Brune if (fas->monitor) { 898794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 899794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 900d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 901eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 902d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 903794bee33SPeter Brune } 904646217ecSPeter Brune for (k = 0; k < fas->max_up_it; k++) { 905646217ecSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 906cc05f883SPeter Brune if (fas->monitor) { 907794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 908794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 909d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 910eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 911d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 912794bee33SPeter Brune } 913646217ecSPeter Brune } 914c90fad12SPeter Brune } else if (snes->pc) { 915c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 9162f7ea302SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 9172f7ea302SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 9182f7ea302SPeter Brune snes->reason = SNES_DIVERGED_INNER; 9192f7ea302SPeter Brune PetscFunctionReturn(0); 9202f7ea302SPeter Brune } 921fe6f9142SPeter Brune } 922fa9694d7SPeter Brune if (fas->next) { 923ee78dd50SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 924794bee33SPeter Brune 925c90fad12SPeter Brune X_c = fas->next->vec_sol; 926293a7e31SPeter Brune Xo_c = fas->next->work[0]; 927c90fad12SPeter Brune F_c = fas->next->vec_func; 928742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 929efe1f98aSPeter Brune 930efe1f98aSPeter Brune /* inject the solution */ 931efe1f98aSPeter Brune if (fas->inject) { 932a3cb90a9SPeter Brune ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 933efe1f98aSPeter Brune } else { 934a3cb90a9SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 935a3cb90a9SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 936efe1f98aSPeter Brune } 937293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 938293a7e31SPeter Brune 939293a7e31SPeter Brune /* restrict the defect */ 940293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 941293a7e31SPeter Brune 942c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 943ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 944c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 945742fe5e2SPeter Brune 946293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 947c90fad12SPeter Brune 948ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 949ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 950c90fad12SPeter Brune 951c90fad12SPeter Brune /* recurse to the next level */ 952f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 953742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 954742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 955742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 956742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 957742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 958742fe5e2SPeter Brune PetscFunctionReturn(0); 959742fe5e2SPeter Brune } 960742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 961ee78dd50SPeter Brune 962fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 963fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 964ee78dd50SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 965293a7e31SPeter Brune } 966742fe5e2SPeter Brune /* up-smooth -- just update using the up-smoother */ 967c90fad12SPeter Brune if (fas->level != 0) { 968d1adcc6fSPeter Brune if (fas->upsmooth) { 969d1adcc6fSPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 970742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 971742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 972742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 973742fe5e2SPeter Brune PetscFunctionReturn(0); 974742fe5e2SPeter Brune } 975d28543b3SPeter Brune } else if (snes->usegs && snes->ops->computegs) { 976794bee33SPeter Brune if (fas->monitor) { 977794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 978794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 979d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 980eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 981d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 982794bee33SPeter Brune } 983646217ecSPeter Brune for (k = 0; k < fas->max_down_it; k++) { 984646217ecSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 985794bee33SPeter Brune if (fas->monitor) { 986794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 987794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 988d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 989eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 990d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 991794bee33SPeter Brune } 992646217ecSPeter Brune } 993c90fad12SPeter Brune } else if (snes->pc) { 994c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 9952f7ea302SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 9962f7ea302SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 9972f7ea302SPeter Brune snes->reason = SNES_DIVERGED_INNER; 9982f7ea302SPeter Brune PetscFunctionReturn(0); 9992f7ea302SPeter Brune } 1000c90fad12SPeter Brune } 1001fe6f9142SPeter Brune } 1002fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1003fa9694d7SPeter Brune 1004fa9694d7SPeter Brune PetscFunctionReturn(0); 1005421d9b32SPeter Brune } 1006421d9b32SPeter Brune 1007421d9b32SPeter Brune #undef __FUNCT__ 1008421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1009421d9b32SPeter Brune 1010421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1011421d9b32SPeter Brune { 1012fa9694d7SPeter Brune PetscErrorCode ierr; 1013fe6f9142SPeter Brune PetscInt i, maxits; 1014f5a6d4f9SBarry Smith Vec X, B,F; 1015fe6f9142SPeter Brune PetscReal fnorm; 1016421d9b32SPeter Brune PetscFunctionBegin; 1017fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1018fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1019fa9694d7SPeter Brune X = snes->vec_sol; 1020fe6f9142SPeter Brune B = snes->vec_rhs; 1021f5a6d4f9SBarry Smith F = snes->vec_func; 1022293a7e31SPeter Brune 1023293a7e31SPeter Brune /*norm setup */ 1024fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1025fe6f9142SPeter Brune snes->iter = 0; 1026fe6f9142SPeter Brune snes->norm = 0.; 1027fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1028fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1029fe6f9142SPeter Brune if (snes->domainerror) { 1030fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1031fe6f9142SPeter Brune PetscFunctionReturn(0); 1032fe6f9142SPeter Brune } 1033fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1034fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1035fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1036fe6f9142SPeter Brune snes->norm = fnorm; 1037fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1038fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1039fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1040fe6f9142SPeter Brune 1041fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1042fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1043fe6f9142SPeter Brune /* test convergence */ 1044fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1045fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1046fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1047fe6f9142SPeter Brune /* Call general purpose update function */ 1048646217ecSPeter Brune 1049fe6f9142SPeter Brune if (snes->ops->update) { 1050fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1051fe6f9142SPeter Brune } 1052646217ecSPeter Brune ierr = FASCycle_Private(snes, X);CHKERRQ(ierr); 1053742fe5e2SPeter Brune 1054742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1055742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1056742fe5e2SPeter Brune PetscFunctionReturn(0); 1057742fe5e2SPeter Brune } 1058c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1059c90fad12SPeter Brune /* Monitor convergence */ 1060c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1061c90fad12SPeter Brune snes->iter = i+1; 1062c90fad12SPeter Brune snes->norm = fnorm; 1063c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1064c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1065c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1066c90fad12SPeter Brune /* Test for convergence */ 1067c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1068c90fad12SPeter Brune if (snes->reason) break; 1069fe6f9142SPeter Brune } 1070fe6f9142SPeter Brune if (i == maxits) { 1071fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1072fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1073fe6f9142SPeter Brune } 1074421d9b32SPeter Brune PetscFunctionReturn(0); 1075421d9b32SPeter Brune } 1076