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); 18421d9b32SPeter Brune 19421d9b32SPeter Brune EXTERN_C_BEGIN 20421d9b32SPeter Brune 21421d9b32SPeter Brune #undef __FUNCT__ 22421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 23421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes) 24421d9b32SPeter Brune { 25421d9b32SPeter Brune SNES_FAS * fas; 26421d9b32SPeter Brune PetscErrorCode ierr; 27421d9b32SPeter Brune 28421d9b32SPeter Brune PetscFunctionBegin; 29421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 30421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 31421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 32421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 33421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 34421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 35421d9b32SPeter Brune 36ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 37ed020824SBarry Smith snes->usespc = PETSC_FALSE; 38ed020824SBarry Smith 39421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 40421d9b32SPeter Brune snes->data = (void*) fas; 41421d9b32SPeter Brune fas->level = 0; 42293a7e31SPeter Brune fas->levels = 1; 43ee78dd50SPeter Brune fas->n_cycles = 1; 44ee78dd50SPeter Brune fas->max_up_it = 1; 45ee78dd50SPeter Brune fas->max_down_it = 1; 46ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 47ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 48421d9b32SPeter Brune fas->next = PETSC_NULL; 49421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 50421d9b32SPeter Brune fas->restrct = PETSC_NULL; 51efe1f98aSPeter Brune fas->inject = PETSC_NULL; 52cc05f883SPeter Brune fas->monitor = PETSC_NULL; 53cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 54eff52c0eSPeter Brune fas->useGS = PETSC_FALSE; 55421d9b32SPeter Brune PetscFunctionReturn(0); 56421d9b32SPeter Brune } 57421d9b32SPeter Brune EXTERN_C_END 58421d9b32SPeter Brune 59421d9b32SPeter Brune #undef __FUNCT__ 60421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 61c79ef259SPeter Brune /*@ 62722262beSPeter Brune SNESFASGetLevels - Gets the number of levels in a FAS. 63c79ef259SPeter Brune 64c79ef259SPeter Brune Input Parameter: 65c79ef259SPeter Brune . snes - the preconditioner context 66c79ef259SPeter Brune 67c79ef259SPeter Brune Output parameter: 68c79ef259SPeter Brune . levels - the number of levels 69c79ef259SPeter Brune 70c79ef259SPeter Brune Level: advanced 71c79ef259SPeter Brune 72c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 73c79ef259SPeter Brune 74c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 75c79ef259SPeter Brune @*/ 76421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 77421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 78421d9b32SPeter Brune PetscFunctionBegin; 79ee1fd11aSPeter Brune *levels = fas->levels; 80421d9b32SPeter Brune PetscFunctionReturn(0); 81421d9b32SPeter Brune } 82421d9b32SPeter Brune 83421d9b32SPeter Brune #undef __FUNCT__ 84646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 85c79ef259SPeter Brune /*@ 86c79ef259SPeter Brune SNESFASSetCycles - Sets the type cycles to use. Use SNESFASSetCyclesOnLevel() for more 87c79ef259SPeter Brune complicated cycling. 88c79ef259SPeter Brune 89c79ef259SPeter Brune Logically Collective on SNES 90c79ef259SPeter Brune 91c79ef259SPeter Brune Input Parameters: 92c79ef259SPeter Brune + snes - the multigrid context 93c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 94c79ef259SPeter Brune 95c79ef259SPeter Brune Options Database Key: 96c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 97c79ef259SPeter Brune 98c79ef259SPeter Brune Level: advanced 99c79ef259SPeter Brune 100c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 101c79ef259SPeter Brune 102c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 103c79ef259SPeter Brune @*/ 104646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 105646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 106646217ecSPeter Brune PetscErrorCode ierr; 107646217ecSPeter Brune PetscFunctionBegin; 108646217ecSPeter Brune fas->n_cycles = cycles; 109646217ecSPeter Brune if (fas->next) { 110646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 111646217ecSPeter Brune } 112646217ecSPeter Brune PetscFunctionReturn(0); 113646217ecSPeter Brune } 114646217ecSPeter Brune 115eff52c0eSPeter Brune #undef __FUNCT__ 116c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 117c79ef259SPeter Brune /*@ 118722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 119c79ef259SPeter Brune 120c79ef259SPeter Brune Logically Collective on SNES 121c79ef259SPeter Brune 122c79ef259SPeter Brune Input Parameters: 123c79ef259SPeter Brune + snes - the multigrid context 124c79ef259SPeter Brune . level - the level to set the number of cycles on 125c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 126c79ef259SPeter Brune 127c79ef259SPeter Brune Level: advanced 128c79ef259SPeter Brune 129c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 130c79ef259SPeter Brune 131c79ef259SPeter Brune .seealso: SNESFASSetCycles() 132c79ef259SPeter Brune @*/ 133c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 134c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 135c79ef259SPeter Brune PetscInt top_level = fas->level,i; 136c79ef259SPeter Brune 137c79ef259SPeter Brune PetscFunctionBegin; 138c79ef259SPeter Brune if (level > top_level) 139c79ef259SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 140c79ef259SPeter Brune /* get to the correct level */ 141c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 142c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 143c79ef259SPeter Brune } 144c79ef259SPeter Brune if (fas->level != level) 145c79ef259SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 146c79ef259SPeter Brune fas->n_cycles = cycles; 147c79ef259SPeter Brune PetscFunctionReturn(0); 148c79ef259SPeter Brune } 149c79ef259SPeter Brune 150c79ef259SPeter Brune #undef __FUNCT__ 151eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 152aeed3662SMatthew G Knepley /*@C 153c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 154c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 155c79ef259SPeter Brune and nonlinear preconditioners. 156c79ef259SPeter Brune 157c79ef259SPeter Brune Logically Collective on SNES 158c79ef259SPeter Brune 159c79ef259SPeter Brune Input Parameters: 160c79ef259SPeter Brune + snes - the multigrid context 161c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 162c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 163c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 164c79ef259SPeter Brune 165c79ef259SPeter Brune Level: advanced 166c79ef259SPeter Brune 167c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 168c79ef259SPeter Brune 169c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 170c79ef259SPeter Brune @*/ 171eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 172eff52c0eSPeter Brune PetscErrorCode ierr = 0; 173eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 174eff52c0eSPeter Brune PetscFunctionBegin; 175eff52c0eSPeter Brune 176eff52c0eSPeter Brune /* use or don't use it according to user wishes*/ 177eff52c0eSPeter Brune fas->useGS = use_gs; 178eff52c0eSPeter Brune if (gsfunc) { 179eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 180eff52c0eSPeter Brune /* push the provided GS up the tree */ 181eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 182eff52c0eSPeter Brune } else if (snes->ops->computegs) { 183eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 184eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 185eff52c0eSPeter Brune } else if (use_gs) { 186eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level); 187eff52c0eSPeter Brune } 188eff52c0eSPeter Brune PetscFunctionReturn(0); 189eff52c0eSPeter Brune } 190eff52c0eSPeter Brune 191eff52c0eSPeter Brune #undef __FUNCT__ 192eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 193aeed3662SMatthew G Knepley /*@C 194c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 195c79ef259SPeter Brune 196c79ef259SPeter Brune Logically Collective on SNES 197c79ef259SPeter Brune 198c79ef259SPeter Brune Input Parameters: 199c79ef259SPeter Brune + snes - the multigrid context 200c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 201c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 202c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 203c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 204c79ef259SPeter Brune 205c79ef259SPeter Brune Level: advanced 206c79ef259SPeter Brune 207c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 208c79ef259SPeter Brune 209c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 210c79ef259SPeter Brune @*/ 211eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 212eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 213eff52c0eSPeter Brune PetscErrorCode ierr; 214eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 215eff52c0eSPeter Brune SNES cur_snes = snes; 216eff52c0eSPeter Brune PetscFunctionBegin; 217eff52c0eSPeter Brune if (level > top_level) 218eff52c0eSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level); 219eff52c0eSPeter Brune /* get to the correct level */ 220eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 221eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 222eff52c0eSPeter Brune cur_snes = fas->next; 223eff52c0eSPeter Brune } 224eff52c0eSPeter Brune if (fas->level != level) 225eff52c0eSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel"); 226eff52c0eSPeter Brune fas->useGS = use_gs; 227eff52c0eSPeter Brune if (gsfunc) { 228eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 229eff52c0eSPeter Brune } 230eff52c0eSPeter Brune PetscFunctionReturn(0); 231eff52c0eSPeter Brune } 232eff52c0eSPeter Brune 233646217ecSPeter Brune #undef __FUNCT__ 234421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 235c79ef259SPeter Brune /*@ 236c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 237c79ef259SPeter Brune level of the FAS hierarchy. 238c79ef259SPeter Brune 239c79ef259SPeter Brune Input Parameters: 240c79ef259SPeter Brune + snes - the multigrid context 241c79ef259SPeter Brune level - the level to get 242c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 243c79ef259SPeter Brune 244c79ef259SPeter Brune Level: advanced 245c79ef259SPeter Brune 246c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 247c79ef259SPeter Brune 248c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 249c79ef259SPeter Brune @*/ 250421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) { 251421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 252421d9b32SPeter Brune PetscInt levels = fas->level; 253421d9b32SPeter Brune PetscInt i; 254421d9b32SPeter Brune PetscFunctionBegin; 255421d9b32SPeter Brune *lsnes = snes; 256421d9b32SPeter Brune if (fas->level < level) { 257421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level."); 258421d9b32SPeter Brune } 259421d9b32SPeter Brune if (level > levels - 1) { 260421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS."); 261421d9b32SPeter Brune } 262421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 263421d9b32SPeter Brune *lsnes = fas->next; 264421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 265421d9b32SPeter Brune } 266421d9b32SPeter Brune if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!"); 267421d9b32SPeter Brune PetscFunctionReturn(0); 268421d9b32SPeter Brune } 269421d9b32SPeter Brune 270421d9b32SPeter Brune #undef __FUNCT__ 271421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 272aeed3662SMatthew G Knepley /*@C 273c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 274c79ef259SPeter Brune Must be called before any other FAS routine. 275c79ef259SPeter Brune 276c79ef259SPeter Brune Input Parameters: 277c79ef259SPeter Brune + snes - the snes context 278c79ef259SPeter Brune . levels - the number of levels 279c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 280c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 281c79ef259SPeter Brune Fortran. 282c79ef259SPeter Brune 283c79ef259SPeter Brune Level: intermediate 284c79ef259SPeter Brune 285c79ef259SPeter Brune Notes: 286c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 287c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 288c79ef259SPeter Brune 289c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 290c79ef259SPeter Brune 291c79ef259SPeter Brune .seealso: SNESFASGetLevels() 292c79ef259SPeter Brune @*/ 293421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 294421d9b32SPeter Brune PetscErrorCode ierr; 295421d9b32SPeter Brune PetscInt i; 296421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 297421d9b32SPeter Brune MPI_Comm comm; 298421d9b32SPeter Brune PetscFunctionBegin; 299ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 300ee1fd11aSPeter Brune if (levels == fas->levels) { 301ee1fd11aSPeter Brune if (!comms) 302ee1fd11aSPeter Brune PetscFunctionReturn(0); 303ee1fd11aSPeter Brune } 304421d9b32SPeter Brune /* user has changed the number of levels; reset */ 305421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 306421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 307421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 308ee1fd11aSPeter Brune fas->next = PETSC_NULL; 309421d9b32SPeter Brune /* setup the finest level */ 310421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 311421d9b32SPeter Brune if (comms) comm = comms[i]; 312421d9b32SPeter Brune fas->level = i; 313421d9b32SPeter Brune fas->levels = levels; 314ee1fd11aSPeter Brune fas->next = PETSC_NULL; 315421d9b32SPeter Brune if (i > 0) { 316421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3174a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 318421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 319794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 320421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 321421d9b32SPeter Brune } 322421d9b32SPeter Brune } 323421d9b32SPeter Brune PetscFunctionReturn(0); 324421d9b32SPeter Brune } 325421d9b32SPeter Brune 326421d9b32SPeter Brune #undef __FUNCT__ 327c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 328c79ef259SPeter Brune /*@ 329c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 330c79ef259SPeter Brune use on all levels. 331c79ef259SPeter Brune 332c79ef259SPeter Brune Logically Collective on PC 333c79ef259SPeter Brune 334c79ef259SPeter Brune Input Parameters: 335c79ef259SPeter Brune + snes - the multigrid context 336c79ef259SPeter Brune - n - the number of smoothing steps 337c79ef259SPeter Brune 338c79ef259SPeter Brune Options Database Key: 339c79ef259SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 340c79ef259SPeter Brune 341c79ef259SPeter Brune Level: advanced 342c79ef259SPeter Brune 343c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 344c79ef259SPeter Brune 345c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 346c79ef259SPeter Brune @*/ 347c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 348c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 349c79ef259SPeter Brune PetscErrorCode ierr = 0; 350c79ef259SPeter Brune PetscFunctionBegin; 351c79ef259SPeter Brune fas->max_up_it = n; 352c79ef259SPeter Brune if (fas->next) { 353c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 354c79ef259SPeter Brune } 355c79ef259SPeter Brune PetscFunctionReturn(0); 356c79ef259SPeter Brune } 357c79ef259SPeter Brune 358c79ef259SPeter Brune #undef __FUNCT__ 359c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 360c79ef259SPeter Brune /*@ 361c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 362c79ef259SPeter Brune use on all levels. 363c79ef259SPeter Brune 364c79ef259SPeter Brune Logically Collective on PC 365c79ef259SPeter Brune 366c79ef259SPeter Brune Input Parameters: 367c79ef259SPeter Brune + snes - the multigrid context 368c79ef259SPeter Brune - n - the number of smoothing steps 369c79ef259SPeter Brune 370c79ef259SPeter Brune Options Database Key: 371c79ef259SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 372c79ef259SPeter Brune 373c79ef259SPeter Brune Level: advanced 374c79ef259SPeter Brune 375c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 376c79ef259SPeter Brune 377c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 378c79ef259SPeter Brune @*/ 379c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 380c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 381c79ef259SPeter Brune PetscErrorCode ierr = 0; 382c79ef259SPeter Brune PetscFunctionBegin; 383c79ef259SPeter Brune fas->max_down_it = n; 384c79ef259SPeter Brune if (fas->next) { 385c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 386c79ef259SPeter Brune } 387c79ef259SPeter Brune PetscFunctionReturn(0); 388c79ef259SPeter Brune } 389c79ef259SPeter Brune 390c79ef259SPeter Brune #undef __FUNCT__ 391421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 392c79ef259SPeter Brune /*@ 393c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 394c79ef259SPeter Brune interpolation from l-1 to the lth level 395c79ef259SPeter Brune 396c79ef259SPeter Brune Input Parameters: 397c79ef259SPeter Brune + snes - the multigrid context 398c79ef259SPeter Brune . mat - the interpolation operator 399c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 400c79ef259SPeter Brune 401c79ef259SPeter Brune Level: advanced 402c79ef259SPeter Brune 403c79ef259SPeter Brune Notes: 404c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 405c79ef259SPeter Brune for the same level. 406c79ef259SPeter Brune 407c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 408c79ef259SPeter Brune out from the matrix size which one it is. 409c79ef259SPeter Brune 410c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 411c79ef259SPeter Brune 412c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale() 413c79ef259SPeter Brune @*/ 414421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 415421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 416d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 417421d9b32SPeter Brune 418421d9b32SPeter Brune PetscFunctionBegin; 419421d9b32SPeter Brune if (level > top_level) 420421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level); 421421d9b32SPeter Brune /* get to the correct level */ 422d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 423421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 424421d9b32SPeter Brune } 425421d9b32SPeter Brune if (fas->level != level) 426421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation"); 427421d9b32SPeter Brune fas->interpolate = mat; 428421d9b32SPeter Brune PetscFunctionReturn(0); 429421d9b32SPeter Brune } 430421d9b32SPeter Brune 431421d9b32SPeter Brune #undef __FUNCT__ 432421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 433c79ef259SPeter Brune /*@ 434c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 435c79ef259SPeter Brune from level l to l-1. 436c79ef259SPeter Brune 437c79ef259SPeter Brune Input Parameters: 438c79ef259SPeter Brune + snes - the multigrid context 439c79ef259SPeter Brune . mat - the restriction matrix 440c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 441c79ef259SPeter Brune 442c79ef259SPeter Brune Level: advanced 443c79ef259SPeter Brune 444c79ef259SPeter Brune Notes: 445c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 446c79ef259SPeter Brune for the same level. 447c79ef259SPeter Brune 448c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 449c79ef259SPeter Brune out from the matrix size which one it is. 450c79ef259SPeter Brune 451c79ef259SPeter Brune If you do not set this, the transpose of the Mat set with PCMGSetInterpolation() 452c79ef259SPeter Brune is used. 453c79ef259SPeter Brune 454c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 455c79ef259SPeter Brune 456c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 457c79ef259SPeter Brune @*/ 458421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 459421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 460d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 461421d9b32SPeter Brune 462421d9b32SPeter Brune PetscFunctionBegin; 463421d9b32SPeter Brune if (level > top_level) 464421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 465421d9b32SPeter Brune /* get to the correct level */ 466d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 467421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 468421d9b32SPeter Brune } 469421d9b32SPeter Brune if (fas->level != level) 470421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 471421d9b32SPeter Brune fas->restrct = mat; 472421d9b32SPeter Brune PetscFunctionReturn(0); 473421d9b32SPeter Brune } 474421d9b32SPeter Brune 475efe1f98aSPeter Brune 476421d9b32SPeter Brune #undef __FUNCT__ 477421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 478c79ef259SPeter Brune /*@ 479c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 480c79ef259SPeter Brune operator from level l to l-1. 481c79ef259SPeter Brune 482c79ef259SPeter Brune Input Parameters: 483c79ef259SPeter Brune + snes - the multigrid context 484c79ef259SPeter Brune . rscale - the restriction scaling 485c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 486c79ef259SPeter Brune 487c79ef259SPeter Brune Level: advanced 488c79ef259SPeter Brune 489c79ef259SPeter Brune Notes: 490c79ef259SPeter Brune This is only used in the case that the injection is not set. 491c79ef259SPeter Brune 492c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 493c79ef259SPeter Brune 494c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 495c79ef259SPeter Brune @*/ 496421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 497421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 498d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 499421d9b32SPeter Brune 500421d9b32SPeter Brune PetscFunctionBegin; 501421d9b32SPeter Brune if (level > top_level) 502421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level); 503421d9b32SPeter Brune /* get to the correct level */ 504d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 505421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 506421d9b32SPeter Brune } 507421d9b32SPeter Brune if (fas->level != level) 508421d9b32SPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction"); 509421d9b32SPeter Brune fas->rscale = rscale; 510421d9b32SPeter Brune PetscFunctionReturn(0); 511421d9b32SPeter Brune } 512421d9b32SPeter Brune 513efe1f98aSPeter Brune 514efe1f98aSPeter Brune #undef __FUNCT__ 515efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 516c79ef259SPeter Brune /*@ 517c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 518c79ef259SPeter Brune from level l to l-1. 519c79ef259SPeter Brune 520c79ef259SPeter Brune Input Parameters: 521c79ef259SPeter Brune + snes - the multigrid context 522c79ef259SPeter Brune . mat - the restriction matrix 523c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 524c79ef259SPeter Brune 525c79ef259SPeter Brune Level: advanced 526c79ef259SPeter Brune 527c79ef259SPeter Brune Notes: 528c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 529c79ef259SPeter Brune project the solution instead. 530c79ef259SPeter Brune 531c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 532c79ef259SPeter Brune 533c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 534c79ef259SPeter Brune @*/ 535efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 536efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 537efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 538efe1f98aSPeter Brune 539efe1f98aSPeter Brune PetscFunctionBegin; 540efe1f98aSPeter Brune if (level > top_level) 541efe1f98aSPeter Brune SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level); 542efe1f98aSPeter Brune /* get to the correct level */ 543efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 544efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 545efe1f98aSPeter Brune } 546efe1f98aSPeter Brune if (fas->level != level) 547efe1f98aSPeter Brune SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection"); 548efe1f98aSPeter Brune fas->inject = mat; 549efe1f98aSPeter Brune PetscFunctionReturn(0); 550efe1f98aSPeter Brune } 551efe1f98aSPeter Brune 552421d9b32SPeter Brune #undef __FUNCT__ 553421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 554421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 555421d9b32SPeter Brune { 55677df8cc4SPeter Brune PetscErrorCode ierr = 0; 557421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 558421d9b32SPeter Brune 559421d9b32SPeter Brune PetscFunctionBegin; 560742fe5e2SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 561742fe5e2SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 562742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 563421d9b32SPeter Brune PetscFunctionReturn(0); 564421d9b32SPeter Brune } 565421d9b32SPeter Brune 566421d9b32SPeter Brune #undef __FUNCT__ 567421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 568421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 569421d9b32SPeter Brune { 570421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 571742fe5e2SPeter Brune PetscErrorCode ierr = 0; 572421d9b32SPeter Brune 573421d9b32SPeter Brune PetscFunctionBegin; 574421d9b32SPeter Brune /* recursively resets and then destroys */ 57551e86f29SPeter Brune if (fas->upsmooth) ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 57651e86f29SPeter Brune if (fas->downsmooth) ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 577efe1f98aSPeter Brune if (fas->inject) { 578efe1f98aSPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 579efe1f98aSPeter Brune } 58051e86f29SPeter Brune if (fas->interpolate == fas->restrct) { 58151e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 58251e86f29SPeter Brune fas->restrct = PETSC_NULL; 58351e86f29SPeter Brune } else { 58451e86f29SPeter Brune if (fas->interpolate) ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 58551e86f29SPeter Brune if (fas->restrct) ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 58651e86f29SPeter Brune } 58751e86f29SPeter Brune if (fas->rscale) ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 588421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 589421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 590ee1fd11aSPeter Brune 591421d9b32SPeter Brune PetscFunctionReturn(0); 592421d9b32SPeter Brune } 593421d9b32SPeter Brune 594421d9b32SPeter Brune #undef __FUNCT__ 595421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 596421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 597421d9b32SPeter Brune { 59848bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 599421d9b32SPeter Brune PetscErrorCode ierr; 600efe1f98aSPeter Brune VecScatter injscatter; 601d1adcc6fSPeter Brune PetscInt dm_levels; 602d1adcc6fSPeter Brune 603421d9b32SPeter Brune 604421d9b32SPeter Brune PetscFunctionBegin; 605eff52c0eSPeter Brune 606cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 607d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 608d1adcc6fSPeter Brune dm_levels++; 609cc05f883SPeter Brune if (dm_levels > fas->levels) { 610d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 611cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 612d1adcc6fSPeter Brune } 613d1adcc6fSPeter Brune } 614d1adcc6fSPeter Brune 615742fe5e2SPeter Brune ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 616cc05f883SPeter Brune 61748bfdf8aSPeter Brune /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 61848bfdf8aSPeter Brune if (fas->upsmooth) { 61948bfdf8aSPeter Brune ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 62048bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 62148bfdf8aSPeter Brune ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 62248bfdf8aSPeter Brune } 62348bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 62448bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 62548bfdf8aSPeter Brune } 62648bfdf8aSPeter Brune if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 62748bfdf8aSPeter Brune ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 62848bfdf8aSPeter Brune } 62948bfdf8aSPeter Brune } 63048bfdf8aSPeter Brune if (fas->downsmooth) { 63148bfdf8aSPeter Brune ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 63248bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 63348bfdf8aSPeter Brune ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 63448bfdf8aSPeter Brune } 63548bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 63648bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 63748bfdf8aSPeter Brune } 63848bfdf8aSPeter Brune if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 63948bfdf8aSPeter Brune ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 64048bfdf8aSPeter Brune } 6411a266240SBarry Smith } 64248bfdf8aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 643646217ecSPeter Brune if (fas->next) { 64448bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->next->ops->computefunction) { 64548bfdf8aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 64648bfdf8aSPeter Brune } 64748bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 64848bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 64948bfdf8aSPeter Brune } 65048bfdf8aSPeter Brune if (snes->ops->computegs && !fas->next->ops->computegs) { 651646217ecSPeter Brune ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 652646217ecSPeter Brune } 653646217ecSPeter Brune } 654421d9b32SPeter Brune if (snes->dm) { 655421d9b32SPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 656421d9b32SPeter Brune if (fas->next) { 657421d9b32SPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 658ee78dd50SPeter Brune if (!fas->next->dm) { 659ee78dd50SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 660ee78dd50SPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 661ee78dd50SPeter Brune } 662421d9b32SPeter Brune /* set the interpolation and restriction from the DM */ 663ee78dd50SPeter Brune if (!fas->interpolate) { 664e727c939SJed Brown ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 665421d9b32SPeter Brune fas->restrct = fas->interpolate; 666421d9b32SPeter Brune } 667efe1f98aSPeter Brune /* set the injection from the DM */ 668efe1f98aSPeter Brune if (!fas->inject) { 669e727c939SJed Brown ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 670efe1f98aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 671efe1f98aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 672efe1f98aSPeter Brune } 673421d9b32SPeter Brune } 674ee78dd50SPeter Brune /* set the DMs of the pre and post-smoothers here */ 675*6bed07a3SBarry Smith if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 676*6bed07a3SBarry Smith if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 677421d9b32SPeter Brune } 678cc05f883SPeter Brune 679fa9694d7SPeter Brune if (fas->next) { 680fa9694d7SPeter Brune /* gotta set up the solution vector for this to work */ 681*6bed07a3SBarry Smith if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);} 682742fe5e2SPeter Brune if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);} 683fa9694d7SPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 684fa9694d7SPeter Brune } 685fa9694d7SPeter Brune /* got to set them all up at once */ 686421d9b32SPeter Brune PetscFunctionReturn(0); 687421d9b32SPeter Brune } 688421d9b32SPeter Brune 689421d9b32SPeter Brune #undef __FUNCT__ 690421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 691421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 692421d9b32SPeter Brune { 693ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 694ee78dd50SPeter Brune PetscInt levels = 1; 695d1adcc6fSPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg; 696421d9b32SPeter Brune PetscErrorCode ierr; 697ee78dd50SPeter Brune const char *def_smooth = SNESNRICHARDSON; 698ee78dd50SPeter Brune char pre_type[256]; 699ee78dd50SPeter Brune char post_type[256]; 700ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 701421d9b32SPeter Brune 702421d9b32SPeter Brune PetscFunctionBegin; 703c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 704ee78dd50SPeter Brune 705ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 706ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 707ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 708c732cbdbSBarry Smith if (!flg && snes->dm) { 709c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 710c732cbdbSBarry Smith levels++; 711d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 712c732cbdbSBarry Smith } 713ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 714ee78dd50SPeter Brune } 715ee78dd50SPeter Brune 716ee78dd50SPeter Brune /* type of pre and/or post smoothers -- set both at once */ 717ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr); 718ee78dd50SPeter Brune ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr); 719d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 720d1adcc6fSPeter Brune if (smoothflg) { 721ee78dd50SPeter Brune ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr); 722ee78dd50SPeter Brune } else { 723d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothup_type", "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr); 724d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr); 725ee78dd50SPeter Brune } 726ee78dd50SPeter Brune 727ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 728d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 729d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 730d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 731ee78dd50SPeter Brune 732646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 733ee78dd50SPeter Brune 734eff52c0eSPeter Brune ierr = PetscOptionsBool("-snes_fas_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr); 735eff52c0eSPeter Brune 736ee78dd50SPeter Brune /* other options for the coarsest level */ 737ee78dd50SPeter Brune if (fas->level == 0) { 738d1adcc6fSPeter Brune ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr); 739eff52c0eSPeter Brune ierr = PetscOptionsBool("-snes_fas_coarse_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr); 740ee78dd50SPeter Brune } 741ee78dd50SPeter Brune 742421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 7438cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 744eff52c0eSPeter Brune 745eff52c0eSPeter Brune if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !fas->useGS)) { 7468cc86e31SPeter Brune const char *prefix; 7478cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 7488cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 7498cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 7508cc86e31SPeter Brune if (fas->level || (fas->levels == 1)) { 751eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr); 7528cc86e31SPeter Brune } else { 7538cc86e31SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 7548cc86e31SPeter Brune } 7558cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 7568cc86e31SPeter Brune ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr); 7578cc86e31SPeter Brune } 7588cc86e31SPeter Brune 759eff52c0eSPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !fas->useGS)) { 76067339d5cSBarry Smith const char *prefix; 76167339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 762ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 76367339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 764eff52c0eSPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr); 765293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 766ee78dd50SPeter Brune ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr); 767ee78dd50SPeter Brune } 7681a266240SBarry Smith if (fas->upsmooth) { 7691a266240SBarry Smith ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr); 7701a266240SBarry Smith } 7711a266240SBarry Smith 7721a266240SBarry Smith if (fas->downsmooth) { 7731a266240SBarry Smith ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr); 7741a266240SBarry Smith } 775ee78dd50SPeter Brune 776742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 777742fe5e2SPeter Brune ierr = SNESSetTolerances(snes, 0.0, 0.0, 0.0, fas->n_cycles, 1000);CHKERRQ(ierr); 778742fe5e2SPeter Brune } 779742fe5e2SPeter Brune 780ee78dd50SPeter Brune if (monflg) { 781646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 782794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 783742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 784293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 785293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 786ee78dd50SPeter Brune } 787ee78dd50SPeter Brune 788ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 7891a266240SBarry Smith if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);} 790421d9b32SPeter Brune PetscFunctionReturn(0); 791421d9b32SPeter Brune } 792421d9b32SPeter Brune 793421d9b32SPeter Brune #undef __FUNCT__ 794421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 795421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 796421d9b32SPeter Brune { 797421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 798421d9b32SPeter Brune PetscBool iascii; 799421d9b32SPeter Brune PetscErrorCode ierr; 800421d9b32SPeter Brune PetscInt levels = fas->levels; 801421d9b32SPeter Brune PetscInt i; 802421d9b32SPeter Brune 803421d9b32SPeter Brune PetscFunctionBegin; 804421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 805421d9b32SPeter Brune if (iascii) { 806421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n", fas->levels);CHKERRQ(ierr); 807421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 808421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 809421d9b32SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n", fas->level);CHKERRQ(ierr); 810ee78dd50SPeter Brune if (fas->upsmooth) { 811ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 812421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 813ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 814421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 815421d9b32SPeter Brune } else { 816ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 817421d9b32SPeter Brune } 818ee78dd50SPeter Brune if (fas->downsmooth) { 819ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 820421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 821ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 822421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 823421d9b32SPeter Brune } else { 824ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 825421d9b32SPeter Brune } 826eff52c0eSPeter Brune if (fas->useGS) { 827eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n", fas->level);CHKERRQ(ierr); 828eff52c0eSPeter Brune } 829421d9b32SPeter Brune if (fas->next) fas = (SNES_FAS *)fas->next->data; 830421d9b32SPeter Brune } 831421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 832421d9b32SPeter Brune } else { 833421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 834421d9b32SPeter Brune } 835421d9b32SPeter Brune PetscFunctionReturn(0); 836421d9b32SPeter Brune } 837421d9b32SPeter Brune 838421d9b32SPeter Brune /* 839fa9694d7SPeter Brune 840fa9694d7SPeter Brune Defines the FAS cycle as: 841fa9694d7SPeter Brune 842fa9694d7SPeter Brune fine problem: F(x) = 0 843fa9694d7SPeter Brune coarse problem: F^c(x) = b^c 844fa9694d7SPeter Brune 845fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 846fa9694d7SPeter Brune 847fa9694d7SPeter Brune correction: 848fa9694d7SPeter Brune 849fa9694d7SPeter Brune x = x + I(x^c - Rx) 850fa9694d7SPeter Brune 851421d9b32SPeter Brune */ 852fa9694d7SPeter Brune 853421d9b32SPeter Brune #undef __FUNCT__ 854421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private" 855646217ecSPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec X) { 856fa9694d7SPeter Brune 857fa9694d7SPeter Brune PetscErrorCode ierr; 858646217ecSPeter Brune Vec X_c, Xo_c, F_c, B_c,F,B; 859fa9694d7SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 860742fe5e2SPeter Brune PetscInt k; 8612d15c84fSPeter Brune PetscReal fnorm; 862742fe5e2SPeter Brune SNESConvergedReason reason; 863421d9b32SPeter Brune 864421d9b32SPeter Brune PetscFunctionBegin; 865f5a6d4f9SBarry Smith F = snes->vec_func; 866646217ecSPeter Brune B = snes->vec_rhs; 867fa9694d7SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 868d1adcc6fSPeter Brune if (fas->downsmooth) { 869d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 870742fe5e2SPeter Brune /* check convergence reason for the smoother */ 871742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 872742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 873742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 874742fe5e2SPeter Brune PetscFunctionReturn(0); 875742fe5e2SPeter Brune } 876d1adcc6fSPeter Brune } else if (snes->ops->computegs) { 877794bee33SPeter Brune if (fas->monitor) { 878794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 879794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 880d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 881eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 882d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 883794bee33SPeter Brune } 884646217ecSPeter Brune for (k = 0; k < fas->max_up_it; k++) { 885646217ecSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 886cc05f883SPeter Brune if (fas->monitor) { 887794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 888794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 889d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 890eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 891d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 892794bee33SPeter Brune } 893646217ecSPeter Brune } 894c90fad12SPeter Brune } else if (snes->pc) { 895c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 896fe6f9142SPeter Brune } 897fa9694d7SPeter Brune if (fas->next) { 898ee78dd50SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 899794bee33SPeter Brune 900c90fad12SPeter Brune X_c = fas->next->vec_sol; 901293a7e31SPeter Brune Xo_c = fas->next->work[0]; 902c90fad12SPeter Brune F_c = fas->next->vec_func; 903742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 904efe1f98aSPeter Brune 905efe1f98aSPeter Brune /* inject the solution */ 906efe1f98aSPeter Brune if (fas->inject) { 907a3cb90a9SPeter Brune ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr); 908efe1f98aSPeter Brune } else { 909a3cb90a9SPeter Brune ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr); 910a3cb90a9SPeter Brune ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr); 911efe1f98aSPeter Brune } 912293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 913293a7e31SPeter Brune 914293a7e31SPeter Brune /* restrict the defect */ 915293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 916293a7e31SPeter Brune 917c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 918ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 919c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 920742fe5e2SPeter Brune 921293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 922c90fad12SPeter Brune 923ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 924ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 925c90fad12SPeter Brune 926c90fad12SPeter Brune /* recurse to the next level */ 927f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 928742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 929742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 930742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 931742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 932742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 933742fe5e2SPeter Brune PetscFunctionReturn(0); 934742fe5e2SPeter Brune } 935742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 936ee78dd50SPeter Brune 937fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 938fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 939ee78dd50SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr); 940293a7e31SPeter Brune } 941742fe5e2SPeter Brune /* up-smooth -- just update using the up-smoother */ 942c90fad12SPeter Brune if (fas->level != 0) { 943d1adcc6fSPeter Brune if (fas->upsmooth) { 944d1adcc6fSPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 945742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 946742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 947742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 948742fe5e2SPeter Brune PetscFunctionReturn(0); 949742fe5e2SPeter Brune } 950d1adcc6fSPeter Brune } else if (snes->ops->computegs) { 951794bee33SPeter Brune if (fas->monitor) { 952794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 953794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 954d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 955eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr); 956d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 957794bee33SPeter Brune } 958646217ecSPeter Brune for (k = 0; k < fas->max_down_it; k++) { 959646217ecSPeter Brune ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr); 960794bee33SPeter Brune if (fas->monitor) { 961794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 962794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 963d1adcc6fSPeter Brune ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 964eff52c0eSPeter Brune ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr); 965d1adcc6fSPeter Brune ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr); 966794bee33SPeter Brune } 967646217ecSPeter Brune } 968c90fad12SPeter Brune } else if (snes->pc) { 969c90fad12SPeter Brune ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr); 970c90fad12SPeter Brune } 971fe6f9142SPeter Brune } 972fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 973fa9694d7SPeter Brune 974fa9694d7SPeter Brune PetscFunctionReturn(0); 975421d9b32SPeter Brune } 976421d9b32SPeter Brune 977421d9b32SPeter Brune #undef __FUNCT__ 978421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 979421d9b32SPeter Brune 980421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 981421d9b32SPeter Brune { 982fa9694d7SPeter Brune PetscErrorCode ierr; 983fe6f9142SPeter Brune PetscInt i, maxits; 984f5a6d4f9SBarry Smith Vec X, B,F; 985fe6f9142SPeter Brune PetscReal fnorm; 986421d9b32SPeter Brune PetscFunctionBegin; 987fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 988fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 989fa9694d7SPeter Brune X = snes->vec_sol; 990fe6f9142SPeter Brune B = snes->vec_rhs; 991f5a6d4f9SBarry Smith F = snes->vec_func; 992293a7e31SPeter Brune 993293a7e31SPeter Brune /*norm setup */ 994fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 995fe6f9142SPeter Brune snes->iter = 0; 996fe6f9142SPeter Brune snes->norm = 0.; 997fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 998fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 999fe6f9142SPeter Brune if (snes->domainerror) { 1000fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1001fe6f9142SPeter Brune PetscFunctionReturn(0); 1002fe6f9142SPeter Brune } 1003fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1004fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1005fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1006fe6f9142SPeter Brune snes->norm = fnorm; 1007fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1008fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1009fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1010fe6f9142SPeter Brune 1011fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1012fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1013fe6f9142SPeter Brune /* test convergence */ 1014fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1015fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1016fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1017fe6f9142SPeter Brune /* Call general purpose update function */ 1018646217ecSPeter Brune 1019fe6f9142SPeter Brune if (snes->ops->update) { 1020fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1021fe6f9142SPeter Brune } 1022646217ecSPeter Brune ierr = FASCycle_Private(snes, X);CHKERRQ(ierr); 1023742fe5e2SPeter Brune 1024742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1025742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1026742fe5e2SPeter Brune PetscFunctionReturn(0); 1027742fe5e2SPeter Brune } 1028c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1029c90fad12SPeter Brune /* Monitor convergence */ 1030c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1031c90fad12SPeter Brune snes->iter = i+1; 1032c90fad12SPeter Brune snes->norm = fnorm; 1033c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1034c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1035c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1036c90fad12SPeter Brune /* Test for convergence */ 1037c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1038c90fad12SPeter Brune if (snes->reason) break; 1039fe6f9142SPeter Brune } 1040fe6f9142SPeter Brune if (i == maxits) { 1041fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1042fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1043fe6f9142SPeter Brune } 1044421d9b32SPeter Brune PetscFunctionReturn(0); 1045421d9b32SPeter Brune } 1046