1421d9b32SPeter Brune /* Defines the basic SNES object */ 26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3421d9b32SPeter Brune 407144faaSPeter Brune const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0}; 507144faaSPeter Brune 6421d9b32SPeter Brune /*MC 7ef536925SPeter Brune 8ef536925SPeter Brune SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 9421d9b32SPeter Brune 10421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections 11421d9b32SPeter Brune 12421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 13421d9b32SPeter Brune M*/ 14421d9b32SPeter Brune 15421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes); 16421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 17421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 18421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 19421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 20421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 216273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 22421d9b32SPeter Brune 23421d9b32SPeter Brune EXTERN_C_BEGIN 24ddebd997SPeter Brune #undef __FUNCT__ 25ddebd997SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_FAS" 26ddebd997SPeter Brune PetscErrorCode SNESLineSearchSetType_FAS(SNES snes, SNESLineSearchType type) 27ddebd997SPeter Brune { 28ddebd997SPeter Brune PetscErrorCode ierr; 29ddebd997SPeter Brune PetscFunctionBegin; 30ddebd997SPeter Brune 31ddebd997SPeter Brune switch (type) { 32ddebd997SPeter Brune case SNES_LS_BASIC: 33ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr); 34ddebd997SPeter Brune break; 35ddebd997SPeter Brune case SNES_LS_BASIC_NONORMS: 36ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr); 37ddebd997SPeter Brune break; 38ddebd997SPeter Brune case SNES_LS_QUADRATIC: 39ddebd997SPeter Brune ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr); 40ddebd997SPeter Brune break; 41ddebd997SPeter Brune default: 42ddebd997SPeter Brune SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type."); 43ddebd997SPeter Brune break; 44ddebd997SPeter Brune } 45ddebd997SPeter Brune snes->ls_type = type; 46ddebd997SPeter Brune PetscFunctionReturn(0); 47ddebd997SPeter Brune } 48ddebd997SPeter Brune EXTERN_C_END 49ddebd997SPeter Brune 50ddebd997SPeter Brune EXTERN_C_BEGIN 51421d9b32SPeter Brune 52421d9b32SPeter Brune #undef __FUNCT__ 53421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 54421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes) 55421d9b32SPeter Brune { 56421d9b32SPeter Brune SNES_FAS * fas; 57421d9b32SPeter Brune PetscErrorCode ierr; 58421d9b32SPeter Brune 59421d9b32SPeter Brune PetscFunctionBegin; 60421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 61421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 62421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 63421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 64421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 65421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 66421d9b32SPeter Brune 67ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 68ed020824SBarry Smith snes->usespc = PETSC_FALSE; 69ed020824SBarry Smith 700e444f03SPeter Brune snes->max_funcs = 30000; 710e444f03SPeter Brune snes->max_its = 10000; 720e444f03SPeter Brune 73421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 74421d9b32SPeter Brune snes->data = (void*) fas; 75421d9b32SPeter Brune fas->level = 0; 76293a7e31SPeter Brune fas->levels = 1; 77ee78dd50SPeter Brune fas->n_cycles = 1; 78ee78dd50SPeter Brune fas->max_up_it = 1; 79ee78dd50SPeter Brune fas->max_down_it = 1; 80ee78dd50SPeter Brune fas->upsmooth = PETSC_NULL; 81ee78dd50SPeter Brune fas->downsmooth = PETSC_NULL; 82421d9b32SPeter Brune fas->next = PETSC_NULL; 836273346dSPeter Brune fas->previous = PETSC_NULL; 84421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 85421d9b32SPeter Brune fas->restrct = PETSC_NULL; 86efe1f98aSPeter Brune fas->inject = PETSC_NULL; 87cc05f883SPeter Brune fas->monitor = PETSC_NULL; 88cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 89ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 90ddebd997SPeter Brune 91ddebd997SPeter Brune ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_FAS",SNESLineSearchSetType_FAS);CHKERRQ(ierr); 92ddebd997SPeter Brune ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr); 93ddebd997SPeter Brune 94421d9b32SPeter Brune PetscFunctionReturn(0); 95421d9b32SPeter Brune } 96421d9b32SPeter Brune EXTERN_C_END 97421d9b32SPeter Brune 98421d9b32SPeter Brune #undef __FUNCT__ 99421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels" 100c79ef259SPeter Brune /*@ 1012e8ce248SJed Brown SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids 102c79ef259SPeter Brune 103c79ef259SPeter Brune Input Parameter: 1042e8ce248SJed Brown . snes - the nonlinear solver context 105c79ef259SPeter Brune 106c79ef259SPeter Brune Output parameter: 107c79ef259SPeter Brune . levels - the number of levels 108c79ef259SPeter Brune 109c79ef259SPeter Brune Level: advanced 110c79ef259SPeter Brune 111c79ef259SPeter Brune .keywords: MG, get, levels, multigrid 112c79ef259SPeter Brune 113c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels() 114c79ef259SPeter Brune @*/ 115421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) { 116421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 117421d9b32SPeter Brune PetscFunctionBegin; 118ee1fd11aSPeter Brune *levels = fas->levels; 119421d9b32SPeter Brune PetscFunctionReturn(0); 120421d9b32SPeter Brune } 121421d9b32SPeter Brune 122421d9b32SPeter Brune #undef __FUNCT__ 123646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles" 124c79ef259SPeter Brune /*@ 1252e8ce248SJed Brown SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited. Use SNESFASSetCyclesOnLevel() for more 126c79ef259SPeter Brune complicated cycling. 127c79ef259SPeter Brune 128c79ef259SPeter Brune Logically Collective on SNES 129c79ef259SPeter Brune 130c79ef259SPeter Brune Input Parameters: 131c79ef259SPeter Brune + snes - the multigrid context 132c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 133c79ef259SPeter Brune 134c79ef259SPeter Brune Options Database Key: 135c79ef259SPeter Brune $ -snes_fas_cycles 1 or 2 136c79ef259SPeter Brune 137c79ef259SPeter Brune Level: advanced 138c79ef259SPeter Brune 139c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 140c79ef259SPeter Brune 141c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel() 142c79ef259SPeter Brune @*/ 143646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) { 144646217ecSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 145646217ecSPeter Brune PetscErrorCode ierr; 146646217ecSPeter Brune PetscFunctionBegin; 147646217ecSPeter Brune fas->n_cycles = cycles; 148646217ecSPeter Brune if (fas->next) { 149646217ecSPeter Brune ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr); 150646217ecSPeter Brune } 151646217ecSPeter Brune PetscFunctionReturn(0); 152646217ecSPeter Brune } 153646217ecSPeter Brune 154eff52c0eSPeter Brune #undef __FUNCT__ 155c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel" 156c79ef259SPeter Brune /*@ 157722262beSPeter Brune SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level. 158c79ef259SPeter Brune 159c79ef259SPeter Brune Logically Collective on SNES 160c79ef259SPeter Brune 161c79ef259SPeter Brune Input Parameters: 162c79ef259SPeter Brune + snes - the multigrid context 163c79ef259SPeter Brune . level - the level to set the number of cycles on 164c79ef259SPeter Brune - cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle 165c79ef259SPeter Brune 166c79ef259SPeter Brune Level: advanced 167c79ef259SPeter Brune 168c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid 169c79ef259SPeter Brune 170c79ef259SPeter Brune .seealso: SNESFASSetCycles() 171c79ef259SPeter Brune @*/ 172c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) { 173c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 174c79ef259SPeter Brune PetscInt top_level = fas->level,i; 175c79ef259SPeter Brune 176c79ef259SPeter Brune PetscFunctionBegin; 1772e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D cycle count from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 178c79ef259SPeter Brune /* get to the correct level */ 179c79ef259SPeter Brune for (i = fas->level; i > level; i--) { 180c79ef259SPeter Brune fas = (SNES_FAS *)fas->next->data; 181c79ef259SPeter Brune } 1822e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 183c79ef259SPeter Brune fas->n_cycles = cycles; 184c79ef259SPeter Brune PetscFunctionReturn(0); 185c79ef259SPeter Brune } 186c79ef259SPeter Brune 187c79ef259SPeter Brune #undef __FUNCT__ 188eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS" 189aeed3662SMatthew G Knepley /*@C 190c79ef259SPeter Brune SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used. 191c79ef259SPeter Brune Use SNESFASSetGSOnLevel() for more complicated staging of smoothers 192c79ef259SPeter Brune and nonlinear preconditioners. 193c79ef259SPeter Brune 194c79ef259SPeter Brune Logically Collective on SNES 195c79ef259SPeter Brune 196c79ef259SPeter Brune Input Parameters: 197c79ef259SPeter Brune + snes - the multigrid context 198c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 199c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 200c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 201c79ef259SPeter Brune 202c79ef259SPeter Brune Level: advanced 203c79ef259SPeter Brune 204c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid 205c79ef259SPeter Brune 206c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel() 207c79ef259SPeter Brune @*/ 208eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 209eff52c0eSPeter Brune PetscErrorCode ierr = 0; 210eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 211eff52c0eSPeter Brune PetscFunctionBegin; 212eff52c0eSPeter Brune 213eff52c0eSPeter Brune if (gsfunc) { 214eff52c0eSPeter Brune ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr); 215eff52c0eSPeter Brune /* push the provided GS up the tree */ 216eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr); 217eff52c0eSPeter Brune } else if (snes->ops->computegs) { 218eff52c0eSPeter Brune /* assume that the user has set the GS solver at this level */ 219eff52c0eSPeter Brune if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr); 220eff52c0eSPeter Brune } else if (use_gs) { 2212e8ce248SJed Brown SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %D", fas->level); 222eff52c0eSPeter Brune } 223eff52c0eSPeter Brune PetscFunctionReturn(0); 224eff52c0eSPeter Brune } 225eff52c0eSPeter Brune 226eff52c0eSPeter Brune #undef __FUNCT__ 227eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel" 228aeed3662SMatthew G Knepley /*@C 229c79ef259SPeter Brune SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level. 230c79ef259SPeter Brune 231c79ef259SPeter Brune Logically Collective on SNES 232c79ef259SPeter Brune 233c79ef259SPeter Brune Input Parameters: 234c79ef259SPeter Brune + snes - the multigrid context 235c79ef259SPeter Brune . level - the level to set the nonlinear smoother on 236c79ef259SPeter Brune . gsfunc - the nonlinear smoother function 237c79ef259SPeter Brune . ctx - the user context for the nonlinear smoother 238c79ef259SPeter Brune - use_gs - whether to use the nonlinear smoother or not 239c79ef259SPeter Brune 240c79ef259SPeter Brune Level: advanced 241c79ef259SPeter Brune 242c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 243c79ef259SPeter Brune 244c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS() 245c79ef259SPeter Brune @*/ 246eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) { 247eff52c0eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 248eff52c0eSPeter Brune PetscErrorCode ierr; 249eff52c0eSPeter Brune PetscInt top_level = fas->level,i; 250eff52c0eSPeter Brune SNES cur_snes = snes; 251eff52c0eSPeter Brune PetscFunctionBegin; 2522e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D GS from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 253eff52c0eSPeter Brune /* get to the correct level */ 254eff52c0eSPeter Brune for (i = fas->level; i > level; i--) { 255eff52c0eSPeter Brune fas = (SNES_FAS *)fas->next->data; 256eff52c0eSPeter Brune cur_snes = fas->next; 257eff52c0eSPeter Brune } 2582e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 259eff52c0eSPeter Brune if (gsfunc) { 2606273346dSPeter Brune ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr); 261eff52c0eSPeter Brune } 262eff52c0eSPeter Brune PetscFunctionReturn(0); 263eff52c0eSPeter Brune } 264eff52c0eSPeter Brune 265646217ecSPeter Brune #undef __FUNCT__ 266421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES" 267c79ef259SPeter Brune /*@ 268c79ef259SPeter Brune SNESFASGetSNES - Gets the SNES corresponding to a particular 269c79ef259SPeter Brune level of the FAS hierarchy. 270c79ef259SPeter Brune 271c79ef259SPeter Brune Input Parameters: 272c79ef259SPeter Brune + snes - the multigrid context 273c79ef259SPeter Brune level - the level to get 274c79ef259SPeter Brune - lsnes - whether to use the nonlinear smoother or not 275c79ef259SPeter Brune 276c79ef259SPeter Brune Level: advanced 277c79ef259SPeter Brune 278c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid 279c79ef259SPeter Brune 280c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels() 281c79ef259SPeter Brune @*/ 282421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes,PetscInt level,SNES *lsnes) { 283421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 284421d9b32SPeter Brune PetscInt i; 2852e8ce248SJed Brown 286421d9b32SPeter Brune PetscFunctionBegin; 2872e8ce248SJed Brown if (level > fas->levels-1) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels); 2882e8ce248SJed Brown if (fas->level < level) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS on level %D, must call from SNES at least as fine as level",level,fas->level); 289421d9b32SPeter Brune *lsnes = snes; 290421d9b32SPeter Brune for (i = fas->level; i > level; i--) { 291421d9b32SPeter Brune *lsnes = fas->next; 292421d9b32SPeter Brune fas = (SNES_FAS *)(*lsnes)->data; 293421d9b32SPeter Brune } 2942e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 295421d9b32SPeter Brune PetscFunctionReturn(0); 296421d9b32SPeter Brune } 297421d9b32SPeter Brune 298421d9b32SPeter Brune #undef __FUNCT__ 29907144faaSPeter Brune #define __FUNCT__ "SNESFASSetType" 30007144faaSPeter Brune /*@ 30107144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS. 30207144faaSPeter Brune e 30307144faaSPeter Brune 30407144faaSPeter Brune 30507144faaSPeter Brune @*/ 30607144faaSPeter Brune PetscErrorCode SNESFASSetType(SNES snes,SNESFASType fastype) 30707144faaSPeter Brune { 30807144faaSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 30907144faaSPeter Brune 31007144faaSPeter Brune PetscFunctionBegin; 31107144faaSPeter Brune PetscValidHeaderSpecific(snes,SNES_CLASSID,1); 31207144faaSPeter Brune PetscValidLogicalCollectiveEnum(snes,fastype,2); 31307144faaSPeter Brune fas->fastype = fastype; 31407144faaSPeter Brune PetscFunctionReturn(0); 31507144faaSPeter Brune } 31607144faaSPeter Brune 31707144faaSPeter Brune #undef __FUNCT__ 318421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels" 319aeed3662SMatthew G Knepley /*@C 320c79ef259SPeter Brune SNESFASSetLevels - Sets the number of levels to use with FAS. 321c79ef259SPeter Brune Must be called before any other FAS routine. 322c79ef259SPeter Brune 323c79ef259SPeter Brune Input Parameters: 324c79ef259SPeter Brune + snes - the snes context 325c79ef259SPeter Brune . levels - the number of levels 326c79ef259SPeter Brune - comms - optional communicators for each level; this is to allow solving the coarser 327c79ef259SPeter Brune problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in 328c79ef259SPeter Brune Fortran. 329c79ef259SPeter Brune 330c79ef259SPeter Brune Level: intermediate 331c79ef259SPeter Brune 332c79ef259SPeter Brune Notes: 333c79ef259SPeter Brune If the number of levels is one then the multigrid uses the -fas_levels prefix 334c79ef259SPeter Brune for setting the level options rather than the -fas_coarse prefix. 335c79ef259SPeter Brune 336c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid 337c79ef259SPeter Brune 338c79ef259SPeter Brune .seealso: SNESFASGetLevels() 339c79ef259SPeter Brune @*/ 340421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) { 341421d9b32SPeter Brune PetscErrorCode ierr; 342421d9b32SPeter Brune PetscInt i; 343421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 3446273346dSPeter Brune SNES prevsnes; 345421d9b32SPeter Brune MPI_Comm comm; 346421d9b32SPeter Brune PetscFunctionBegin; 347ee1fd11aSPeter Brune comm = ((PetscObject)snes)->comm; 348ee1fd11aSPeter Brune if (levels == fas->levels) { 349ee1fd11aSPeter Brune if (!comms) 350ee1fd11aSPeter Brune PetscFunctionReturn(0); 351ee1fd11aSPeter Brune } 352421d9b32SPeter Brune /* user has changed the number of levels; reset */ 353421d9b32SPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 354421d9b32SPeter Brune /* destroy any coarser levels if necessary */ 355421d9b32SPeter Brune if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr); 356ee1fd11aSPeter Brune fas->next = PETSC_NULL; 3576273346dSPeter Brune fas->previous = PETSC_NULL; 3586273346dSPeter Brune prevsnes = snes; 359421d9b32SPeter Brune /* setup the finest level */ 360421d9b32SPeter Brune for (i = levels - 1; i >= 0; i--) { 361421d9b32SPeter Brune if (comms) comm = comms[i]; 362421d9b32SPeter Brune fas->level = i; 363421d9b32SPeter Brune fas->levels = levels; 364ee1fd11aSPeter Brune fas->next = PETSC_NULL; 365421d9b32SPeter Brune if (i > 0) { 366421d9b32SPeter Brune ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr); 3674a6b3fb8SBarry Smith ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr); 368421d9b32SPeter Brune ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr); 369794bee33SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr); 3706273346dSPeter Brune ((SNES_FAS *)fas->next->data)->previous = prevsnes; 3716273346dSPeter Brune prevsnes = fas->next; 3726273346dSPeter Brune fas = (SNES_FAS *)prevsnes->data; 373421d9b32SPeter Brune } 374421d9b32SPeter Brune } 375421d9b32SPeter Brune PetscFunctionReturn(0); 376421d9b32SPeter Brune } 377421d9b32SPeter Brune 378421d9b32SPeter Brune #undef __FUNCT__ 379c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp" 380c79ef259SPeter Brune /*@ 381c79ef259SPeter Brune SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to 382c79ef259SPeter Brune use on all levels. 383c79ef259SPeter Brune 384fde0ff24SPeter Brune Logically Collective on SNES 385c79ef259SPeter Brune 386c79ef259SPeter Brune Input Parameters: 387c79ef259SPeter Brune + snes - the multigrid context 388c79ef259SPeter Brune - n - the number of smoothing steps 389c79ef259SPeter Brune 390c79ef259SPeter Brune Options Database Key: 391d28543b3SPeter Brune . -snes_fas_smoothup <n> - Sets number of pre-smoothing steps 392c79ef259SPeter Brune 393c79ef259SPeter Brune Level: advanced 394c79ef259SPeter Brune 395c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 396c79ef259SPeter Brune 397c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown() 398c79ef259SPeter Brune @*/ 399c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) { 400c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 401c79ef259SPeter Brune PetscErrorCode ierr = 0; 402c79ef259SPeter Brune PetscFunctionBegin; 403c79ef259SPeter Brune fas->max_up_it = n; 404c79ef259SPeter Brune if (fas->next) { 405c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr); 406c79ef259SPeter Brune } 407c79ef259SPeter Brune PetscFunctionReturn(0); 408c79ef259SPeter Brune } 409c79ef259SPeter Brune 410c79ef259SPeter Brune #undef __FUNCT__ 411c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown" 412c79ef259SPeter Brune /*@ 413c79ef259SPeter Brune SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to 414c79ef259SPeter Brune use on all levels. 415c79ef259SPeter Brune 416fde0ff24SPeter Brune Logically Collective on SNES 417c79ef259SPeter Brune 418c79ef259SPeter Brune Input Parameters: 419c79ef259SPeter Brune + snes - the multigrid context 420c79ef259SPeter Brune - n - the number of smoothing steps 421c79ef259SPeter Brune 422c79ef259SPeter Brune Options Database Key: 423d28543b3SPeter Brune . -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps 424c79ef259SPeter Brune 425c79ef259SPeter Brune Level: advanced 426c79ef259SPeter Brune 427c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid 428c79ef259SPeter Brune 429c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp() 430c79ef259SPeter Brune @*/ 431c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) { 432c79ef259SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 433c79ef259SPeter Brune PetscErrorCode ierr = 0; 434c79ef259SPeter Brune PetscFunctionBegin; 435c79ef259SPeter Brune fas->max_down_it = n; 436c79ef259SPeter Brune if (fas->next) { 437c79ef259SPeter Brune ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr); 438c79ef259SPeter Brune } 439c79ef259SPeter Brune PetscFunctionReturn(0); 440c79ef259SPeter Brune } 441c79ef259SPeter Brune 442c79ef259SPeter Brune #undef __FUNCT__ 443421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation" 444c79ef259SPeter Brune /*@ 445c79ef259SPeter Brune SNESFASSetInterpolation - Sets the function to be used to calculate the 446c79ef259SPeter Brune interpolation from l-1 to the lth level 447c79ef259SPeter Brune 448c79ef259SPeter Brune Input Parameters: 449c79ef259SPeter Brune + snes - the multigrid context 450c79ef259SPeter Brune . mat - the interpolation operator 451c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [do not supply 0] 452c79ef259SPeter Brune 453c79ef259SPeter Brune Level: advanced 454c79ef259SPeter Brune 455c79ef259SPeter Brune Notes: 456c79ef259SPeter Brune Usually this is the same matrix used also to set the restriction 457c79ef259SPeter Brune for the same level. 458c79ef259SPeter Brune 459c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 460c79ef259SPeter Brune out from the matrix size which one it is. 461c79ef259SPeter Brune 462c79ef259SPeter Brune .keywords: FAS, multigrid, set, interpolate, level 463c79ef259SPeter Brune 464bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale() 465c79ef259SPeter Brune @*/ 466421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) { 467421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 468d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 469bccf9bb3SJed Brown PetscErrorCode ierr; 470421d9b32SPeter Brune 471421d9b32SPeter Brune PetscFunctionBegin; 4722e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D interpolation from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 473421d9b32SPeter Brune /* get to the correct level */ 474d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 475421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 476421d9b32SPeter Brune } 4772e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 478bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 479bd4e12b0SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 480421d9b32SPeter Brune fas->interpolate = mat; 481421d9b32SPeter Brune PetscFunctionReturn(0); 482421d9b32SPeter Brune } 483421d9b32SPeter Brune 484421d9b32SPeter Brune #undef __FUNCT__ 485421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction" 486c79ef259SPeter Brune /*@ 487c79ef259SPeter Brune SNESFASSetRestriction - Sets the function to be used to restrict the defect 488c79ef259SPeter Brune from level l to l-1. 489c79ef259SPeter Brune 490c79ef259SPeter Brune Input Parameters: 491c79ef259SPeter Brune + snes - the multigrid context 492c79ef259SPeter Brune . mat - the restriction matrix 493c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 494c79ef259SPeter Brune 495c79ef259SPeter Brune Level: advanced 496c79ef259SPeter Brune 497c79ef259SPeter Brune Notes: 498c79ef259SPeter Brune Usually this is the same matrix used also to set the interpolation 499c79ef259SPeter Brune for the same level. 500c79ef259SPeter Brune 501c79ef259SPeter Brune One can pass in the interpolation matrix or its transpose; PETSc figures 502c79ef259SPeter Brune out from the matrix size which one it is. 503c79ef259SPeter Brune 504fde0ff24SPeter Brune If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation() 505c79ef259SPeter Brune is used. 506c79ef259SPeter Brune 507c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 508c79ef259SPeter Brune 509c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection() 510c79ef259SPeter Brune @*/ 511421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) { 512421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 513d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 514bccf9bb3SJed Brown PetscErrorCode ierr; 515421d9b32SPeter Brune 516421d9b32SPeter Brune PetscFunctionBegin; 5172e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D restriction from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 518421d9b32SPeter Brune /* get to the correct level */ 519d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 520421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 521421d9b32SPeter Brune } 5222e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 523bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 524bd4e12b0SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 525421d9b32SPeter Brune fas->restrct = mat; 526421d9b32SPeter Brune PetscFunctionReturn(0); 527421d9b32SPeter Brune } 528421d9b32SPeter Brune 529421d9b32SPeter Brune #undef __FUNCT__ 530421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale" 531c79ef259SPeter Brune /*@ 532c79ef259SPeter Brune SNESFASSetRScale - Sets the scaling factor of the restriction 533c79ef259SPeter Brune operator from level l to l-1. 534c79ef259SPeter Brune 535c79ef259SPeter Brune Input Parameters: 536c79ef259SPeter Brune + snes - the multigrid context 537c79ef259SPeter Brune . rscale - the restriction scaling 538c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 539c79ef259SPeter Brune 540c79ef259SPeter Brune Level: advanced 541c79ef259SPeter Brune 542c79ef259SPeter Brune Notes: 543c79ef259SPeter Brune This is only used in the case that the injection is not set. 544c79ef259SPeter Brune 545c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 546c79ef259SPeter Brune 547c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction() 548c79ef259SPeter Brune @*/ 549421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) { 550421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 551d9ce41dcSBarry Smith PetscInt top_level = fas->level,i; 552bccf9bb3SJed Brown PetscErrorCode ierr; 553421d9b32SPeter Brune 554421d9b32SPeter Brune PetscFunctionBegin; 5552e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D rscale from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 556421d9b32SPeter Brune /* get to the correct level */ 557d9ce41dcSBarry Smith for (i = fas->level; i > level; i--) { 558421d9b32SPeter Brune fas = (SNES_FAS *)fas->next->data; 559421d9b32SPeter Brune } 5602e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 561bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr); 562bd4e12b0SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 563421d9b32SPeter Brune fas->rscale = rscale; 564421d9b32SPeter Brune PetscFunctionReturn(0); 565421d9b32SPeter Brune } 566421d9b32SPeter Brune 567efe1f98aSPeter Brune 568efe1f98aSPeter Brune #undef __FUNCT__ 569efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection" 570c79ef259SPeter Brune /*@ 571c79ef259SPeter Brune SNESFASSetInjection - Sets the function to be used to inject the solution 572c79ef259SPeter Brune from level l to l-1. 573c79ef259SPeter Brune 574c79ef259SPeter Brune Input Parameters: 575c79ef259SPeter Brune + snes - the multigrid context 576c79ef259SPeter Brune . mat - the restriction matrix 577c79ef259SPeter Brune - level - the level (0 is coarsest) to supply [Do not supply 0] 578c79ef259SPeter Brune 579c79ef259SPeter Brune Level: advanced 580c79ef259SPeter Brune 581c79ef259SPeter Brune Notes: 582c79ef259SPeter Brune If you do not set this, the restriction and rscale is used to 583c79ef259SPeter Brune project the solution instead. 584c79ef259SPeter Brune 585c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level 586c79ef259SPeter Brune 587c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction() 588c79ef259SPeter Brune @*/ 589efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) { 590efe1f98aSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 591efe1f98aSPeter Brune PetscInt top_level = fas->level,i; 592bccf9bb3SJed Brown PetscErrorCode ierr; 593efe1f98aSPeter Brune 594efe1f98aSPeter Brune PetscFunctionBegin; 5952e8ce248SJed Brown if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D injection from level %D of SNESFAS with %D levels total",level,top_level,fas->levels); 596efe1f98aSPeter Brune /* get to the correct level */ 597efe1f98aSPeter Brune for (i = fas->level; i > level; i--) { 598efe1f98aSPeter Brune fas = (SNES_FAS *)fas->next->data; 599efe1f98aSPeter Brune } 6002e8ce248SJed Brown if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt"); 601bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr); 602bd4e12b0SJed Brown ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 603efe1f98aSPeter Brune fas->inject = mat; 604efe1f98aSPeter Brune PetscFunctionReturn(0); 605efe1f98aSPeter Brune } 606efe1f98aSPeter Brune 607421d9b32SPeter Brune #undef __FUNCT__ 608421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 609421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 610421d9b32SPeter Brune { 61177df8cc4SPeter Brune PetscErrorCode ierr = 0; 612421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 613421d9b32SPeter Brune 614421d9b32SPeter Brune PetscFunctionBegin; 615bccf9bb3SJed Brown ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr); 616bccf9bb3SJed Brown ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr); 6173dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 618bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 619bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 620bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 6213dccd265SPeter Brune 6223dccd265SPeter Brune if (fas->upsmooth) ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr); 6233dccd265SPeter Brune if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr); 624742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 62509dc8347SJed Brown ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","",PETSC_NULL);CHKERRQ(ierr); 626421d9b32SPeter Brune PetscFunctionReturn(0); 627421d9b32SPeter Brune } 628421d9b32SPeter Brune 629421d9b32SPeter Brune #undef __FUNCT__ 630421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 631421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 632421d9b32SPeter Brune { 633421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 634742fe5e2SPeter Brune PetscErrorCode ierr = 0; 635421d9b32SPeter Brune 636421d9b32SPeter Brune PetscFunctionBegin; 637421d9b32SPeter Brune /* recursively resets and then destroys */ 63879d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 639421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 640421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 641ee1fd11aSPeter Brune 642421d9b32SPeter Brune PetscFunctionReturn(0); 643421d9b32SPeter Brune } 644421d9b32SPeter Brune 645421d9b32SPeter Brune #undef __FUNCT__ 646421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 647421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 648421d9b32SPeter Brune { 64948bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 650421d9b32SPeter Brune PetscErrorCode ierr; 651efe1f98aSPeter Brune VecScatter injscatter; 652d1adcc6fSPeter Brune PetscInt dm_levels; 6533dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 654d1adcc6fSPeter Brune 655421d9b32SPeter Brune PetscFunctionBegin; 656eff52c0eSPeter Brune 657cc05f883SPeter Brune if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) { 658d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 659d1adcc6fSPeter Brune dm_levels++; 660cc05f883SPeter Brune if (dm_levels > fas->levels) { 6613dccd265SPeter Brune 6622e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 6633dccd265SPeter Brune vec_sol = snes->vec_sol; 6643dccd265SPeter Brune vec_func = snes->vec_func; 6653dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 6663dccd265SPeter Brune vec_rhs = snes->vec_rhs; 6673dccd265SPeter Brune snes->vec_sol = PETSC_NULL; 6683dccd265SPeter Brune snes->vec_func = PETSC_NULL; 6693dccd265SPeter Brune snes->vec_sol_update = PETSC_NULL; 6703dccd265SPeter Brune snes->vec_rhs = PETSC_NULL; 6713dccd265SPeter Brune 6723dccd265SPeter Brune /* reset the number of levels */ 673d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 674cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 6753dccd265SPeter Brune 6763dccd265SPeter Brune snes->vec_sol = vec_sol; 6773dccd265SPeter Brune snes->vec_func = vec_func; 6783dccd265SPeter Brune snes->vec_rhs = vec_rhs; 6793dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 680d1adcc6fSPeter Brune } 681d1adcc6fSPeter Brune } 682d1adcc6fSPeter Brune 6833dccd265SPeter Brune if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 6843dccd265SPeter Brune 68507144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 686fde0ff24SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 68707144faaSPeter Brune } else { 688ddebd997SPeter Brune ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 68907144faaSPeter Brune } 690cc05f883SPeter Brune 69179d9a41aSPeter Brune if (snes->dm) { 69279d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 69379d9a41aSPeter Brune if (fas->next) { 69479d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 69579d9a41aSPeter Brune if (!fas->next->dm) { 69679d9a41aSPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr); 69779d9a41aSPeter Brune ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr); 69879d9a41aSPeter Brune } 69979d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 70079d9a41aSPeter Brune if (!fas->interpolate) { 70179d9a41aSPeter Brune ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 702bccf9bb3SJed Brown if (!fas->restrct) { 703bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 70479d9a41aSPeter Brune fas->restrct = fas->interpolate; 70579d9a41aSPeter Brune } 706bccf9bb3SJed Brown } 70779d9a41aSPeter Brune /* set the injection from the DM */ 70879d9a41aSPeter Brune if (!fas->inject) { 70979d9a41aSPeter Brune ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 71079d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 71179d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 71279d9a41aSPeter Brune } 71379d9a41aSPeter Brune } 71479d9a41aSPeter Brune /* set the DMs of the pre and post-smoothers here */ 71579d9a41aSPeter Brune if (fas->upsmooth) {ierr = SNESSetDM(fas->upsmooth, snes->dm);CHKERRQ(ierr);} 71679d9a41aSPeter Brune if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);} 71779d9a41aSPeter Brune } 71879d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 71979d9a41aSPeter Brune 72079d9a41aSPeter Brune if (fas->next) { 72179d9a41aSPeter Brune if (fas->galerkin) { 72279d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr); 72379d9a41aSPeter Brune } else { 72479d9a41aSPeter Brune if (snes->ops->computefunction && !fas->next->ops->computefunction) { 72579d9a41aSPeter Brune ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 72679d9a41aSPeter Brune } 72779d9a41aSPeter Brune if (snes->ops->computejacobian && !fas->next->ops->computejacobian) { 72879d9a41aSPeter Brune ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 72979d9a41aSPeter Brune } 73079d9a41aSPeter Brune if (snes->ops->computegs && !fas->next->ops->computegs) { 73179d9a41aSPeter Brune ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 73279d9a41aSPeter Brune } 73379d9a41aSPeter Brune } 73479d9a41aSPeter Brune } 73579d9a41aSPeter Brune 73679d9a41aSPeter Brune if (fas->next) { 73779d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 738938e4a01SJed Brown if (!fas->next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_sol);CHKERRQ(ierr);} 739938e4a01SJed Brown if (!fas->next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_rhs);CHKERRQ(ierr);} 74079d9a41aSPeter Brune ierr = SNESSetUp(fas->next);CHKERRQ(ierr); 74179d9a41aSPeter Brune } 74279d9a41aSPeter Brune 74348bfdf8aSPeter Brune /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */ 74448bfdf8aSPeter Brune if (fas->upsmooth) { 74548bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) { 74648bfdf8aSPeter Brune ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 74748bfdf8aSPeter Brune } 74848bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) { 74948bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 75048bfdf8aSPeter Brune } 75148bfdf8aSPeter Brune if (snes->ops->computegs && !fas->upsmooth->ops->computegs) { 75248bfdf8aSPeter Brune ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 75348bfdf8aSPeter Brune } 754fde0ff24SPeter Brune ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr); 75548bfdf8aSPeter Brune } 75648bfdf8aSPeter Brune if (fas->downsmooth) { 75748bfdf8aSPeter Brune if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) { 75848bfdf8aSPeter Brune ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr); 75948bfdf8aSPeter Brune } 76048bfdf8aSPeter Brune if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) { 76148bfdf8aSPeter Brune ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr); 76248bfdf8aSPeter Brune } 76348bfdf8aSPeter Brune if (snes->ops->computegs && !fas->downsmooth->ops->computegs) { 76448bfdf8aSPeter Brune ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr); 76548bfdf8aSPeter Brune } 766fde0ff24SPeter Brune ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr); 7671a266240SBarry Smith } 768cc05f883SPeter Brune 7696273346dSPeter Brune /* setup FAS work vectors */ 7706273346dSPeter Brune if (fas->galerkin) { 7716273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 7726273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 7736273346dSPeter Brune } 7746273346dSPeter Brune 77579d9a41aSPeter Brune 776fa9694d7SPeter Brune /* got to set them all up at once */ 777421d9b32SPeter Brune PetscFunctionReturn(0); 778421d9b32SPeter Brune } 779421d9b32SPeter Brune 780421d9b32SPeter Brune #undef __FUNCT__ 781421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 782421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 783421d9b32SPeter Brune { 784ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 785ee78dd50SPeter Brune PetscInt levels = 1; 786fde0ff24SPeter Brune PetscBool flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg; 787421d9b32SPeter Brune PetscErrorCode ierr; 788ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 78907144faaSPeter Brune SNESFASType fastype; 790fde0ff24SPeter Brune const char *optionsprefix; 791fde0ff24SPeter Brune const char *prefix; 792421d9b32SPeter Brune 793421d9b32SPeter Brune PetscFunctionBegin; 794c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 795ee78dd50SPeter Brune 796ee78dd50SPeter Brune /* number of levels -- only process on the finest level */ 797ee78dd50SPeter Brune if (fas->levels - 1 == fas->level) { 798ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 799c732cbdbSBarry Smith if (!flg && snes->dm) { 800c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 801c732cbdbSBarry Smith levels++; 802d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 803c732cbdbSBarry Smith } 804ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 805ee78dd50SPeter Brune } 80607144faaSPeter Brune fastype = fas->fastype; 80707144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 80807144faaSPeter Brune if (flg) { 80907144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 81007144faaSPeter Brune } 811ee78dd50SPeter Brune 812fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 813fde0ff24SPeter Brune 814fde0ff24SPeter Brune /* smoother setup options */ 815fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr); 816fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr); 817fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr); 818fde0ff24SPeter Brune if (fas->level == 0) { 819fde0ff24SPeter Brune ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr); 820fde0ff24SPeter Brune } 821ee78dd50SPeter Brune /* options for the number of preconditioning cycles and cycle type */ 822fde0ff24SPeter Brune 823d1adcc6fSPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr); 824ee78dd50SPeter Brune 8256273346dSPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr); 8266273346dSPeter Brune 827646217ecSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 828ee78dd50SPeter Brune 829ee78dd50SPeter Brune 830421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 8318cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 832eff52c0eSPeter Brune 833fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothcoarseflg) { 834fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 835fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 836fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 837fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr); 838fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 839fde0ff24SPeter Brune } 840fde0ff24SPeter Brune 841fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothdownflg) { 8428cc86e31SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 8438cc86e31SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 8448cc86e31SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 845fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr); 8468cc86e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 8478cc86e31SPeter Brune } 8488cc86e31SPeter Brune 849fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) { 85067339d5cSBarry Smith ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 851ee78dd50SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 85267339d5cSBarry Smith ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 853fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr); 854293a7e31SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 855ee78dd50SPeter Brune } 856fde0ff24SPeter Brune 857fde0ff24SPeter Brune if ((!fas->downsmooth) && smoothflg) { 858fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 859fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr); 860fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr); 861fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr); 862fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 863fde0ff24SPeter Brune } 864fde0ff24SPeter Brune 865fde0ff24SPeter Brune if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) { 866fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr); 867fde0ff24SPeter Brune ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr); 868fde0ff24SPeter Brune ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr); 869fde0ff24SPeter Brune ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr); 870fde0ff24SPeter Brune ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr); 871fde0ff24SPeter Brune } 872fde0ff24SPeter Brune 8731a266240SBarry Smith if (fas->upsmooth) { 874fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8751a266240SBarry Smith } 8761a266240SBarry Smith 8771a266240SBarry Smith if (fas->downsmooth) { 878fde0ff24SPeter Brune ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr); 8791a266240SBarry Smith } 880ee78dd50SPeter Brune 881742fe5e2SPeter Brune if (fas->level != fas->levels - 1) { 882fde0ff24SPeter Brune ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr); 883742fe5e2SPeter Brune } 884742fe5e2SPeter Brune 885ce11c761SPeter Brune /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */ 886ce11c761SPeter Brune if (!fas->downsmooth) { 88793dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 888*5d115551SPeter Brune ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 889ce11c761SPeter Brune if (fas->level == 0) { 89079d9a41aSPeter Brune ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr); 891ce11c761SPeter Brune } 892ce11c761SPeter Brune } 893ce11c761SPeter Brune 894ce11c761SPeter Brune if (!fas->upsmooth) { 89593dfaa4eSPeter Brune ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 896*5d115551SPeter Brune ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr); 897ce11c761SPeter Brune } 898ce11c761SPeter Brune 899ee78dd50SPeter Brune if (monflg) { 900646217ecSPeter Brune fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr); 901794bee33SPeter Brune /* set the monitors for the upsmoother and downsmoother */ 9022f7ea302SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 903742fe5e2SPeter Brune ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 904293a7e31SPeter Brune if (fas->upsmooth) ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 905293a7e31SPeter Brune if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr); 906d28543b3SPeter Brune } else { 907d28543b3SPeter Brune /* unset the monitors on the coarse levels */ 908d28543b3SPeter Brune if (fas->level != fas->levels - 1) { 909d28543b3SPeter Brune ierr = SNESMonitorCancel(snes);CHKERRQ(ierr); 910d28543b3SPeter Brune } 911ee78dd50SPeter Brune } 912ee78dd50SPeter Brune 913ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 914d28543b3SPeter Brune if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);} 915421d9b32SPeter Brune PetscFunctionReturn(0); 916421d9b32SPeter Brune } 917421d9b32SPeter Brune 918421d9b32SPeter Brune #undef __FUNCT__ 919421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 920421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 921421d9b32SPeter Brune { 922421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 923421d9b32SPeter Brune PetscBool iascii; 924421d9b32SPeter Brune PetscErrorCode ierr; 925421d9b32SPeter Brune 926421d9b32SPeter Brune PetscFunctionBegin; 927421d9b32SPeter Brune ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr); 928421d9b32SPeter Brune if (iascii) { 9292e8ce248SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n", fas->levels);CHKERRQ(ierr); 930421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 931bd4e12b0SJed Brown ierr = PetscViewerASCIIPrintf(viewer, "level: %D\n", fas->level);CHKERRQ(ierr); 932ee78dd50SPeter Brune if (fas->upsmooth) { 93339a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 934421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 935ee78dd50SPeter Brune ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr); 936421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 937421d9b32SPeter Brune } else { 938ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n", fas->level);CHKERRQ(ierr); 939421d9b32SPeter Brune } 940ee78dd50SPeter Brune if (fas->downsmooth) { 94139a4b097SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n", fas->level);CHKERRQ(ierr); 942421d9b32SPeter Brune ierr = PetscViewerASCIIPushTab(viewer); 943ee78dd50SPeter Brune ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr); 944421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 945421d9b32SPeter Brune } else { 946ee78dd50SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n", fas->level);CHKERRQ(ierr); 947421d9b32SPeter Brune } 948421d9b32SPeter Brune ierr = PetscViewerASCIIPopTab(viewer); 949421d9b32SPeter Brune } else { 950421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 951421d9b32SPeter Brune } 952421d9b32SPeter Brune PetscFunctionReturn(0); 953421d9b32SPeter Brune } 954421d9b32SPeter Brune 955421d9b32SPeter Brune #undef __FUNCT__ 95639bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth" 95739bd7f45SPeter Brune /* 95839bd7f45SPeter Brune Defines the action of the downsmoother 95939bd7f45SPeter Brune */ 96039bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){ 96139bd7f45SPeter Brune PetscErrorCode ierr = 0; 962fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 963742fe5e2SPeter Brune SNESConvergedReason reason; 96439bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 9654b32a720SPeter Brune Vec G, W, Y, FPC; 966fde0ff24SPeter Brune PetscBool lssuccess; 967fde0ff24SPeter Brune PetscInt k; 968421d9b32SPeter Brune PetscFunctionBegin; 969fde0ff24SPeter Brune G = snes->work[1]; 970fde0ff24SPeter Brune W = snes->work[2]; 971fde0ff24SPeter Brune Y = snes->work[3]; 972d1adcc6fSPeter Brune if (fas->downsmooth) { 973d1adcc6fSPeter Brune ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr); 974742fe5e2SPeter Brune /* check convergence reason for the smoother */ 975742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr); 976742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 977742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 978742fe5e2SPeter Brune PetscFunctionReturn(0); 979742fe5e2SPeter Brune } 9804b32a720SPeter Brune ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 9814b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 982fde0ff24SPeter Brune } else { 983fde0ff24SPeter Brune /* basic richardson smoother */ 984fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 985794bee33SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 986794bee33SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 987fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 98805b53524SPeter Brune ierr = SNESLineSearchPreCheckApply(snes,X,Y,PETSC_NULL);CHKERRQ(ierr); 98905b53524SPeter Brune ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 990fde0ff24SPeter Brune if (!lssuccess) { 991fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 992fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 9932f7ea302SPeter Brune PetscFunctionReturn(0); 9942f7ea302SPeter Brune } 995fe6f9142SPeter Brune } 996fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 997fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 998fde0ff24SPeter Brune fnorm = gnorm; 999fde0ff24SPeter Brune } 1000fde0ff24SPeter Brune } 100139bd7f45SPeter Brune PetscFunctionReturn(0); 100239bd7f45SPeter Brune } 100339bd7f45SPeter Brune 100439bd7f45SPeter Brune 100539bd7f45SPeter Brune #undef __FUNCT__ 100639bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth" 100739bd7f45SPeter Brune /* 100807144faaSPeter Brune Defines the action of the upsmoother 100939bd7f45SPeter Brune */ 101039bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) { 101139bd7f45SPeter Brune PetscErrorCode ierr = 0; 1012fde0ff24SPeter Brune PetscReal fnorm, gnorm, ynorm, xnorm = 0.0; 101339bd7f45SPeter Brune SNESConvergedReason reason; 101439bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 10154b32a720SPeter Brune Vec G, W, Y, FPC; 1016fde0ff24SPeter Brune PetscBool lssuccess; 1017fde0ff24SPeter Brune PetscInt k; 101839bd7f45SPeter Brune PetscFunctionBegin; 1019fde0ff24SPeter Brune G = snes->work[1]; 1020fde0ff24SPeter Brune W = snes->work[2]; 1021fde0ff24SPeter Brune Y = snes->work[3]; 102239bd7f45SPeter Brune if (fas->upsmooth) { 1023fde0ff24SPeter Brune ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr); 102439bd7f45SPeter Brune /* check convergence reason for the smoother */ 1025fde0ff24SPeter Brune ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr); 102639bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 102739bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 102839bd7f45SPeter Brune PetscFunctionReturn(0); 102939bd7f45SPeter Brune } 10304b32a720SPeter Brune ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 10314b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 1032fde0ff24SPeter Brune } else { 1033fde0ff24SPeter Brune /* basic richardson smoother */ 1034fde0ff24SPeter Brune for (k = 0; k < fas->max_up_it; k++) { 103539bd7f45SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 103639bd7f45SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1037fde0ff24SPeter Brune ierr = VecCopy(F, Y);CHKERRQ(ierr); 103805b53524SPeter Brune ierr = SNESLineSearchPreCheckApply(snes,X,Y,PETSC_NULL);CHKERRQ(ierr); 103905b53524SPeter Brune ierr = SNESLineSearchApply(snes,X,F,Y,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr); 1040fde0ff24SPeter Brune if (!lssuccess) { 1041fde0ff24SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 1042fde0ff24SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 104339bd7f45SPeter Brune PetscFunctionReturn(0); 104439bd7f45SPeter Brune } 104539bd7f45SPeter Brune } 1046fde0ff24SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1047fde0ff24SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1048fde0ff24SPeter Brune fnorm = gnorm; 1049fde0ff24SPeter Brune } 1050fde0ff24SPeter Brune } 105139bd7f45SPeter Brune PetscFunctionReturn(0); 105239bd7f45SPeter Brune } 105339bd7f45SPeter Brune 105439bd7f45SPeter Brune #undef __FUNCT__ 1055938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 1056938e4a01SJed Brown /*@ 1057938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 1058938e4a01SJed Brown 1059938e4a01SJed Brown Collective 1060938e4a01SJed Brown 1061938e4a01SJed Brown Input Arguments: 1062938e4a01SJed Brown . snes - SNESFAS 1063938e4a01SJed Brown 1064938e4a01SJed Brown Output Arguments: 1065938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 1066938e4a01SJed Brown 1067938e4a01SJed Brown Level: developer 1068938e4a01SJed Brown 1069938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 1070938e4a01SJed Brown @*/ 1071938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 1072938e4a01SJed Brown { 1073938e4a01SJed Brown PetscErrorCode ierr; 1074938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 1075938e4a01SJed Brown 1076938e4a01SJed Brown PetscFunctionBegin; 1077938e4a01SJed Brown if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 1078938e4a01SJed Brown else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 1079938e4a01SJed Brown else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 1080938e4a01SJed Brown PetscFunctionReturn(0); 1081938e4a01SJed Brown } 1082938e4a01SJed Brown 1083e9923e8dSJed Brown #undef __FUNCT__ 1084e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 1085e9923e8dSJed Brown /*@ 1086e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 1087e9923e8dSJed Brown 1088e9923e8dSJed Brown Collective 1089e9923e8dSJed Brown 1090e9923e8dSJed Brown Input Arguments: 1091e9923e8dSJed Brown + fine - SNES from which to restrict 1092e9923e8dSJed Brown - Xfine - vector to restrict 1093e9923e8dSJed Brown 1094e9923e8dSJed Brown Output Arguments: 1095e9923e8dSJed Brown . Xcoarse - result of restriction 1096e9923e8dSJed Brown 1097e9923e8dSJed Brown Level: developer 1098e9923e8dSJed Brown 1099e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 1100e9923e8dSJed Brown @*/ 1101e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 1102e9923e8dSJed Brown { 1103e9923e8dSJed Brown PetscErrorCode ierr; 1104e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 1105e9923e8dSJed Brown 1106e9923e8dSJed Brown PetscFunctionBegin; 1107e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 1108e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 1109e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 1110e9923e8dSJed Brown if (fas->inject) { 1111e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 1112e9923e8dSJed Brown } else { 1113e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 1114e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 1115e9923e8dSJed Brown } 1116e9923e8dSJed Brown PetscFunctionReturn(0); 1117e9923e8dSJed Brown } 1118e9923e8dSJed Brown 1119e9923e8dSJed Brown #undef __FUNCT__ 112039bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection" 112139bd7f45SPeter Brune /* 112239bd7f45SPeter Brune 112339bd7f45SPeter Brune Performs the FAS coarse correction as: 112439bd7f45SPeter Brune 112539bd7f45SPeter Brune fine problem: F(x) = 0 112639bd7f45SPeter Brune coarse problem: F^c(x) = b^c 112739bd7f45SPeter Brune 112839bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 112939bd7f45SPeter Brune 113039bd7f45SPeter Brune with correction: 113139bd7f45SPeter Brune 113239bd7f45SPeter Brune 113339bd7f45SPeter Brune 113439bd7f45SPeter Brune */ 113539a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 113639bd7f45SPeter Brune PetscErrorCode ierr; 113739bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 113839bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 113939bd7f45SPeter Brune SNESConvergedReason reason; 114039bd7f45SPeter Brune PetscFunctionBegin; 1141fa9694d7SPeter Brune if (fas->next) { 1142c90fad12SPeter Brune X_c = fas->next->vec_sol; 1143293a7e31SPeter Brune Xo_c = fas->next->work[0]; 1144c90fad12SPeter Brune F_c = fas->next->vec_func; 1145742fe5e2SPeter Brune B_c = fas->next->vec_rhs; 1146efe1f98aSPeter Brune 1147938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 1148293a7e31SPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 1149293a7e31SPeter Brune 1150293a7e31SPeter Brune /* restrict the defect */ 1151293a7e31SPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 1152293a7e31SPeter Brune 1153c90fad12SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 1154ee78dd50SPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 1155c90fad12SPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 1156742fe5e2SPeter Brune 1157293a7e31SPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 1158c90fad12SPeter Brune 1159ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 1160ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 1161c90fad12SPeter Brune 1162c90fad12SPeter Brune /* recurse to the next level */ 1163f5a6d4f9SBarry Smith fas->next->vec_rhs = B_c; 1164742fe5e2SPeter Brune /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */ 1165742fe5e2SPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 1166742fe5e2SPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 1167742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 1168742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 1169742fe5e2SPeter Brune PetscFunctionReturn(0); 1170742fe5e2SPeter Brune } 1171742fe5e2SPeter Brune /* fas->next->vec_rhs = PETSC_NULL; */ 1172ee78dd50SPeter Brune 1173fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 1174fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 117539bd7f45SPeter Brune ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr); 1176293a7e31SPeter Brune } 117739bd7f45SPeter Brune PetscFunctionReturn(0); 117839bd7f45SPeter Brune } 117939bd7f45SPeter Brune 118039bd7f45SPeter Brune #undef __FUNCT__ 118139bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive" 118239bd7f45SPeter Brune /* 118339bd7f45SPeter Brune 118439bd7f45SPeter Brune The additive cycle looks like: 118539bd7f45SPeter Brune 118607144faaSPeter Brune xhat = x 118707144faaSPeter Brune xhat = dS(x, b) 118807144faaSPeter Brune x = coarsecorrection(xhat, b_d) 118907144faaSPeter Brune x = x + nu*(xhat - x); 119039bd7f45SPeter Brune (optional) x = uS(x, b) 119139bd7f45SPeter Brune 119239bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 119339bd7f45SPeter Brune 119439bd7f45SPeter Brune */ 119539bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) { 119607144faaSPeter Brune Vec F, B, Xhat; 1197ddebd997SPeter Brune Vec X_c, Xo_c, F_c, B_c, G, W; 119839bd7f45SPeter Brune PetscErrorCode ierr; 119907144faaSPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 120007144faaSPeter Brune SNESConvergedReason reason; 1201ddebd997SPeter Brune PetscReal xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.; 1202ddebd997SPeter Brune PetscBool lssucceed; 120339bd7f45SPeter Brune PetscFunctionBegin; 120439bd7f45SPeter Brune 120539bd7f45SPeter Brune F = snes->vec_func; 120639bd7f45SPeter Brune B = snes->vec_rhs; 1207fde0ff24SPeter Brune Xhat = snes->work[3]; 1208fde0ff24SPeter Brune G = snes->work[1]; 1209fde0ff24SPeter Brune W = snes->work[2]; 121007144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 121107144faaSPeter Brune /* recurse first */ 121207144faaSPeter Brune if (fas->next) { 121307144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 121439bd7f45SPeter Brune 121507144faaSPeter Brune X_c = fas->next->vec_sol; 121607144faaSPeter Brune Xo_c = fas->next->work[0]; 121707144faaSPeter Brune F_c = fas->next->vec_func; 121807144faaSPeter Brune B_c = fas->next->vec_rhs; 121939bd7f45SPeter Brune 1220938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 122107144faaSPeter Brune ierr = VecScale(F, -1.0);CHKERRQ(ierr); 122207144faaSPeter Brune 122307144faaSPeter Brune /* restrict the defect */ 122407144faaSPeter Brune ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr); 122507144faaSPeter Brune 122607144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 122707144faaSPeter Brune fas->next->vec_rhs = PETSC_NULL; /*unset the RHS to evaluate function instead of residual*/ 122807144faaSPeter Brune ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr); 122907144faaSPeter Brune 123007144faaSPeter Brune ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr); /* add F_c(X) to the RHS */ 123107144faaSPeter Brune 123207144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 123307144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 123407144faaSPeter Brune 123507144faaSPeter Brune /* recurse */ 123607144faaSPeter Brune fas->next->vec_rhs = B_c; 123707144faaSPeter Brune ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr); 123807144faaSPeter Brune 123907144faaSPeter Brune /* smooth on this level */ 124007144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 124107144faaSPeter Brune 124207144faaSPeter Brune ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr); 124307144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 124407144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 124507144faaSPeter Brune PetscFunctionReturn(0); 124607144faaSPeter Brune } 124707144faaSPeter Brune 124807144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 124907144faaSPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 1250ddebd997SPeter Brune ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr); 125107144faaSPeter Brune 1252ddebd997SPeter Brune /* additive correction of the coarse direction*/ 1253ddebd997SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1254ddebd997SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 1255eb1825c3SPeter Brune ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr); 125605b53524SPeter Brune ierr = SNESLineSearchPreCheckApply(snes, X,Xhat,PETSC_NULL);CHKERRQ(ierr); 125705b53524SPeter Brune ierr = SNESLineSearchApply(snes,X,F,Xhat,fnorm,xnorm,W,G,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr); 1258ddebd997SPeter Brune ierr = VecCopy(W, X);CHKERRQ(ierr); 1259ddebd997SPeter Brune ierr = VecCopy(G, F);CHKERRQ(ierr); 1260ddebd997SPeter Brune fnorm = gnorm; 126107144faaSPeter Brune } else { 126207144faaSPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 126307144faaSPeter Brune } 126439bd7f45SPeter Brune PetscFunctionReturn(0); 126539bd7f45SPeter Brune } 126639bd7f45SPeter Brune 126739bd7f45SPeter Brune #undef __FUNCT__ 126839bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative" 126939bd7f45SPeter Brune /* 127039bd7f45SPeter Brune 127139bd7f45SPeter Brune Defines the FAS cycle as: 127239bd7f45SPeter Brune 127339bd7f45SPeter Brune fine problem: F(x) = 0 127439bd7f45SPeter Brune coarse problem: F^c(x) = b^c 127539bd7f45SPeter Brune 127639bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 127739bd7f45SPeter Brune 127839bd7f45SPeter Brune correction: 127939bd7f45SPeter Brune 128039bd7f45SPeter Brune x = x + I(x^c - Rx) 128139bd7f45SPeter Brune 128239bd7f45SPeter Brune */ 128339bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) { 128439bd7f45SPeter Brune 128539bd7f45SPeter Brune PetscErrorCode ierr; 128639bd7f45SPeter Brune Vec F,B; 128739bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 128839bd7f45SPeter Brune 128939bd7f45SPeter Brune PetscFunctionBegin; 129039bd7f45SPeter Brune F = snes->vec_func; 129139bd7f45SPeter Brune B = snes->vec_rhs; 129239bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 129339bd7f45SPeter Brune ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr); 129439bd7f45SPeter Brune 129539a4b097SPeter Brune ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 129639bd7f45SPeter Brune 1297c90fad12SPeter Brune if (fas->level != 0) { 129839bd7f45SPeter Brune ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr); 1299fe6f9142SPeter Brune } 1300fe6f9142SPeter Brune ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr); 1301fa9694d7SPeter Brune 1302fa9694d7SPeter Brune PetscFunctionReturn(0); 1303421d9b32SPeter Brune } 1304421d9b32SPeter Brune 1305421d9b32SPeter Brune #undef __FUNCT__ 1306421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 1307421d9b32SPeter Brune 1308421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 1309421d9b32SPeter Brune { 1310fa9694d7SPeter Brune PetscErrorCode ierr; 1311fe6f9142SPeter Brune PetscInt i, maxits; 1312ddb5aff1SPeter Brune Vec X, F; 1313fe6f9142SPeter Brune PetscReal fnorm; 1314b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 1315b17ce1afSJed Brown DM dm; 1316b17ce1afSJed Brown 1317421d9b32SPeter Brune PetscFunctionBegin; 1318fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 1319fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 1320fa9694d7SPeter Brune X = snes->vec_sol; 1321f5a6d4f9SBarry Smith F = snes->vec_func; 1322293a7e31SPeter Brune 1323293a7e31SPeter Brune /*norm setup */ 1324fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1325fe6f9142SPeter Brune snes->iter = 0; 1326fe6f9142SPeter Brune snes->norm = 0.; 1327fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1328fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 1329fe6f9142SPeter Brune if (snes->domainerror) { 1330fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 1331fe6f9142SPeter Brune PetscFunctionReturn(0); 1332fe6f9142SPeter Brune } 1333fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1334fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 1335fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1336fe6f9142SPeter Brune snes->norm = fnorm; 1337fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1338fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 1339fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 1340fe6f9142SPeter Brune 1341fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 1342fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 1343fe6f9142SPeter Brune /* test convergence */ 1344fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1345fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 1346b17ce1afSJed Brown 1347b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 1348b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 1349b17ce1afSJed Brown DM dmcoarse; 1350b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 1351b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 1352b17ce1afSJed Brown dm = dmcoarse; 1353b17ce1afSJed Brown } 1354b17ce1afSJed Brown 1355fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 1356fe6f9142SPeter Brune /* Call general purpose update function */ 1357646217ecSPeter Brune 1358fe6f9142SPeter Brune if (snes->ops->update) { 1359fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 1360fe6f9142SPeter Brune } 136107144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 136239bd7f45SPeter Brune ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 136307144faaSPeter Brune } else { 136407144faaSPeter Brune ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr); 136507144faaSPeter Brune } 1366742fe5e2SPeter Brune 1367742fe5e2SPeter Brune /* check for FAS cycle divergence */ 1368742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 1369742fe5e2SPeter Brune PetscFunctionReturn(0); 1370742fe5e2SPeter Brune } 1371c90fad12SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 1372c90fad12SPeter Brune /* Monitor convergence */ 1373c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 1374c90fad12SPeter Brune snes->iter = i+1; 1375c90fad12SPeter Brune snes->norm = fnorm; 1376c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 1377c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 1378c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1379c90fad12SPeter Brune /* Test for convergence */ 1380c90fad12SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1381c90fad12SPeter Brune if (snes->reason) break; 1382fe6f9142SPeter Brune } 1383fe6f9142SPeter Brune if (i == maxits) { 1384fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1385fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1386fe6f9142SPeter Brune } 1387421d9b32SPeter Brune PetscFunctionReturn(0); 1388421d9b32SPeter Brune } 1389