1421d9b32SPeter Brune /* Defines the basic SNES object */ 26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3421d9b32SPeter Brune 46a6fc655SJed Brown const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0}; 507144faaSPeter Brune 6421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes); 7421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 8421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes); 9421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 10421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 11421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 126273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *); 13ab8d36c9SPeter Brune extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *); 14421d9b32SPeter Brune 151fbfccc6SJed Brown /*MC 161fbfccc6SJed Brown 171fbfccc6SJed Brown SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 181fbfccc6SJed Brown 19d3bc2379SPeter Brune The nonlinear problem is solved by correction using coarse versions 20d3bc2379SPeter Brune of the nonlinear problem. This problem is perturbed so that a projected 21d3bc2379SPeter Brune solution of the fine problem elicits no correction from the coarse problem. 22d3bc2379SPeter Brune 23d3bc2379SPeter Brune Options Database: 24d3bc2379SPeter Brune + -snes_fas_levels - The number of levels 25d3bc2379SPeter Brune . -snes_fas_cycles<1> - The number of cycles -- 1 for V, 2 for W 26d3bc2379SPeter Brune . -snes_fas_type<additive, multiplicative> - Additive or multiplicative cycle 27d3bc2379SPeter Brune . -snes_fas_galerkin<PETSC_FALSE> - Form coarse problems by projection back upon the fine problem 28d3bc2379SPeter Brune . -snes_fas_smoothup<1> - The number of iterations of the post-smoother 29d3bc2379SPeter Brune . -snes_fas_smoothdown<1> - The number of iterations of the pre-smoother 30d3bc2379SPeter Brune . -snes_fas_monitor - Monitor progress of all of the levels 31d3bc2379SPeter Brune . -fas_levels_snes_ - SNES options for all smoothers 327d84e935SPeter Brune . -fas_levels_cycle_snes_ - SNES options for all cycles 33d3bc2379SPeter Brune . -fas_levels_i_snes_ - SNES options for the smoothers on level i 347d84e935SPeter Brune . -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i 35d3bc2379SPeter Brune - -fas_coarse_snes_ - SNES options for the coarsest smoother 36d3bc2379SPeter Brune 37d3bc2379SPeter Brune Notes: 38d3bc2379SPeter Brune The organization of the FAS solver is slightly different from the organization of PCMG 39d3bc2379SPeter Brune As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance. 40d3bc2379SPeter Brune The cycle SNES instance may be used for monitoring convergence on a particular level. 411fbfccc6SJed Brown 427d84e935SPeter Brune Level: beginner 431fbfccc6SJed Brown 44d3bc2379SPeter Brune .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 451fbfccc6SJed Brown M*/ 46421d9b32SPeter Brune 47421d9b32SPeter Brune #undef __FUNCT__ 48421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 491fbfccc6SJed Brown PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes) 50421d9b32SPeter Brune { 51421d9b32SPeter Brune SNES_FAS * fas; 52421d9b32SPeter Brune PetscErrorCode ierr; 53421d9b32SPeter Brune 54421d9b32SPeter Brune PetscFunctionBegin; 55421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 56421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 57421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 58421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 59421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 60421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 61421d9b32SPeter Brune 62ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 63ed020824SBarry Smith snes->usespc = PETSC_FALSE; 64ed020824SBarry Smith 6588976e71SPeter Brune if (!snes->tolerancesset) { 660e444f03SPeter Brune snes->max_funcs = 30000; 670e444f03SPeter Brune snes->max_its = 10000; 6888976e71SPeter Brune } 690e444f03SPeter Brune 70421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 71421d9b32SPeter Brune snes->data = (void*) fas; 72421d9b32SPeter Brune fas->level = 0; 73293a7e31SPeter Brune fas->levels = 1; 74ee78dd50SPeter Brune fas->n_cycles = 1; 75ee78dd50SPeter Brune fas->max_up_it = 1; 76ee78dd50SPeter Brune fas->max_down_it = 1; 77ab8d36c9SPeter Brune fas->smoothu = PETSC_NULL; 78ab8d36c9SPeter Brune fas->smoothd = PETSC_NULL; 79421d9b32SPeter Brune fas->next = PETSC_NULL; 806273346dSPeter Brune fas->previous = PETSC_NULL; 81ab8d36c9SPeter Brune fas->fine = snes; 82421d9b32SPeter Brune fas->interpolate = PETSC_NULL; 83421d9b32SPeter Brune fas->restrct = PETSC_NULL; 84efe1f98aSPeter Brune fas->inject = PETSC_NULL; 85cc05f883SPeter Brune fas->monitor = PETSC_NULL; 86cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 87ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 88efe1f98aSPeter Brune PetscFunctionReturn(0); 89efe1f98aSPeter Brune } 90efe1f98aSPeter Brune 91421d9b32SPeter Brune #undef __FUNCT__ 92421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 93421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 94421d9b32SPeter Brune { 9577df8cc4SPeter Brune PetscErrorCode ierr = 0; 96421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 97421d9b32SPeter Brune 98421d9b32SPeter Brune PetscFunctionBegin; 99ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr); 100ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 1013dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 102bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 103bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 104bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 105742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 10622c1e704SPeter Brune 107421d9b32SPeter Brune PetscFunctionReturn(0); 108421d9b32SPeter Brune } 109421d9b32SPeter Brune 110421d9b32SPeter Brune #undef __FUNCT__ 111421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 112421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 113421d9b32SPeter Brune { 114421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS *)snes->data; 115742fe5e2SPeter Brune PetscErrorCode ierr = 0; 116421d9b32SPeter Brune 117421d9b32SPeter Brune PetscFunctionBegin; 118421d9b32SPeter Brune /* recursively resets and then destroys */ 11979d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 120421d9b32SPeter Brune if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 121421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 122ee1fd11aSPeter Brune 123421d9b32SPeter Brune PetscFunctionReturn(0); 124421d9b32SPeter Brune } 125421d9b32SPeter Brune 126421d9b32SPeter Brune #undef __FUNCT__ 127421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 128421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 129421d9b32SPeter Brune { 13048bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 131421d9b32SPeter Brune PetscErrorCode ierr; 132efe1f98aSPeter Brune VecScatter injscatter; 133d1adcc6fSPeter Brune PetscInt dm_levels; 1343dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 135ab8d36c9SPeter Brune SNES next; 136ab8d36c9SPeter Brune PetscBool isFine; 137f89ba88eSPeter Brune SNESLineSearch linesearch; 138f89ba88eSPeter Brune SNESLineSearch slinesearch; 139f89ba88eSPeter Brune void *lsprectx,*lspostctx; 140f89ba88eSPeter Brune SNESLineSearchPreCheckFunc lsprefunc; 141f89ba88eSPeter Brune SNESLineSearchPostCheckFunc lspostfunc; 142421d9b32SPeter Brune PetscFunctionBegin; 143eff52c0eSPeter Brune 144ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 145ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 146d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 147d1adcc6fSPeter Brune dm_levels++; 148cc05f883SPeter Brune if (dm_levels > fas->levels) { 1492e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 1503dccd265SPeter Brune vec_sol = snes->vec_sol; 1513dccd265SPeter Brune vec_func = snes->vec_func; 1523dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 1533dccd265SPeter Brune vec_rhs = snes->vec_rhs; 1543dccd265SPeter Brune snes->vec_sol = PETSC_NULL; 1553dccd265SPeter Brune snes->vec_func = PETSC_NULL; 1563dccd265SPeter Brune snes->vec_sol_update = PETSC_NULL; 1573dccd265SPeter Brune snes->vec_rhs = PETSC_NULL; 1583dccd265SPeter Brune 1593dccd265SPeter Brune /* reset the number of levels */ 160d1adcc6fSPeter Brune ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr); 161cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 1623dccd265SPeter Brune 1633dccd265SPeter Brune snes->vec_sol = vec_sol; 1643dccd265SPeter Brune snes->vec_func = vec_func; 1653dccd265SPeter Brune snes->vec_rhs = vec_rhs; 1663dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 167d1adcc6fSPeter Brune } 168d1adcc6fSPeter Brune } 169ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 170ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 1713dccd265SPeter Brune 17207144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 173e4ed7901SPeter Brune ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 17407144faaSPeter Brune } else { 175e7f468e7SPeter Brune ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 17607144faaSPeter Brune } 177cc05f883SPeter Brune 178ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 179ab8d36c9SPeter Brune if (!fas->smoothd) { 180ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 181ab8d36c9SPeter Brune } 182534ebe21SPeter Brune if (!fas->smoothu && fas->level != fas->levels - 1) { 183534ebe21SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr); 184534ebe21SPeter Brune } 185ab8d36c9SPeter Brune 18679d9a41aSPeter Brune if (snes->dm) { 187ab8d36c9SPeter Brune /* set the smoother DMs properly */ 188ab8d36c9SPeter Brune if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 189ab8d36c9SPeter Brune ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 19079d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 191ab8d36c9SPeter Brune if (next) { 19279d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 193ab8d36c9SPeter Brune if (!next->dm) { 194ab8d36c9SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr); 195ab8d36c9SPeter Brune ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 19679d9a41aSPeter Brune } 19779d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 19879d9a41aSPeter Brune if (!fas->interpolate) { 199ab8d36c9SPeter Brune ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 200bccf9bb3SJed Brown if (!fas->restrct) { 201bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 20279d9a41aSPeter Brune fas->restrct = fas->interpolate; 20379d9a41aSPeter Brune } 204bccf9bb3SJed Brown } 20579d9a41aSPeter Brune /* set the injection from the DM */ 20679d9a41aSPeter Brune if (!fas->inject) { 207ab8d36c9SPeter Brune ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 20879d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 20979d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 21079d9a41aSPeter Brune } 21179d9a41aSPeter Brune } 21279d9a41aSPeter Brune } 21379d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 21479d9a41aSPeter Brune if (fas->galerkin) { 215ab8d36c9SPeter Brune if (next) 216ab8d36c9SPeter Brune ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 217ab8d36c9SPeter Brune if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 218ab8d36c9SPeter Brune if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 21979d9a41aSPeter Brune } 22079d9a41aSPeter Brune 221534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 222534ebe21SPeter Brune if (fas->smoothd){ 223534ebe21SPeter Brune if (fas->level == 0) { 224534ebe21SPeter Brune ierr = SNESSetNormType(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); 225534ebe21SPeter Brune } else { 226534ebe21SPeter Brune ierr = SNESSetNormType(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 227534ebe21SPeter Brune } 2287fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr); 229534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); 230f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 231f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr); 232f89ba88eSPeter Brune ierr = SNESLineSearchGetPreCheck(linesearch,&lsprefunc,&lsprectx);CHKERRQ(ierr); 233f89ba88eSPeter Brune ierr = SNESLineSearchGetPostCheck(linesearch,&lspostfunc,&lspostctx);CHKERRQ(ierr); 234f89ba88eSPeter Brune ierr = SNESLineSearchSetPreCheck(slinesearch,lsprefunc,lsprectx);CHKERRQ(ierr); 235f89ba88eSPeter Brune ierr = SNESLineSearchSetPostCheck(slinesearch,lspostfunc,lspostctx);CHKERRQ(ierr); 236f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 237534ebe21SPeter Brune } 238534ebe21SPeter Brune 239534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 240534ebe21SPeter Brune if (fas->smoothu){ 241534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 242534ebe21SPeter Brune ierr = SNESSetNormType(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); 243534ebe21SPeter Brune } else { 244534ebe21SPeter Brune ierr = SNESSetNormType(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 245534ebe21SPeter Brune } 2467fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr); 247534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); 248f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 249f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr); 250f89ba88eSPeter Brune ierr = SNESLineSearchGetPreCheck(linesearch,&lsprefunc,&lsprectx);CHKERRQ(ierr); 251f89ba88eSPeter Brune ierr = SNESLineSearchGetPostCheck(linesearch,&lspostfunc,&lspostctx);CHKERRQ(ierr); 252f89ba88eSPeter Brune ierr = SNESLineSearchSetPreCheck(slinesearch,lsprefunc,lsprectx);CHKERRQ(ierr); 253f89ba88eSPeter Brune ierr = SNESLineSearchSetPostCheck(slinesearch,lspostfunc,lspostctx);CHKERRQ(ierr); 254f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 255534ebe21SPeter Brune } 256d06165b7SPeter Brune 257ab8d36c9SPeter Brune if (next) { 25879d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 259ab8d36c9SPeter Brune if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 260ab8d36c9SPeter Brune if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 2617fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); 262f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 263f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(fas->next,&slinesearch);CHKERRQ(ierr); 264f89ba88eSPeter Brune ierr = SNESLineSearchGetPreCheck(linesearch,&lsprefunc,&lsprectx);CHKERRQ(ierr); 265f89ba88eSPeter Brune ierr = SNESLineSearchGetPostCheck(linesearch,&lspostfunc,&lspostctx);CHKERRQ(ierr); 266f89ba88eSPeter Brune ierr = SNESLineSearchSetPreCheck(slinesearch,lsprefunc,lsprectx);CHKERRQ(ierr); 267f89ba88eSPeter Brune ierr = SNESLineSearchSetPostCheck(slinesearch,lspostfunc,lspostctx);CHKERRQ(ierr); 268f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 269ab8d36c9SPeter Brune ierr = SNESSetUp(next);CHKERRQ(ierr); 27079d9a41aSPeter Brune } 2716273346dSPeter Brune /* setup FAS work vectors */ 2726273346dSPeter Brune if (fas->galerkin) { 2736273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 2746273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 2756273346dSPeter Brune } 276421d9b32SPeter Brune PetscFunctionReturn(0); 277421d9b32SPeter Brune } 278421d9b32SPeter Brune 279421d9b32SPeter Brune #undef __FUNCT__ 280421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 281421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 282421d9b32SPeter Brune { 283ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 284ee78dd50SPeter Brune PetscInt levels = 1; 2854d26bfa5SPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE; 286421d9b32SPeter Brune PetscErrorCode ierr; 287ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 28807144faaSPeter Brune SNESFASType fastype; 289fde0ff24SPeter Brune const char *optionsprefix; 290f1c6b773SPeter Brune SNESLineSearch linesearch; 29166585501SPeter Brune PetscInt m, n_up, n_down; 292ab8d36c9SPeter Brune SNES next; 293ab8d36c9SPeter Brune PetscBool isFine; 294421d9b32SPeter Brune 295421d9b32SPeter Brune PetscFunctionBegin; 296ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 297c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 298ee78dd50SPeter Brune 299ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 300ab8d36c9SPeter Brune if (isFine) { 301ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 302c732cbdbSBarry Smith if (!flg && snes->dm) { 303c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 304c732cbdbSBarry Smith levels++; 305d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 306c732cbdbSBarry Smith } 307ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 30807144faaSPeter Brune fastype = fas->fastype; 30907144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 31007144faaSPeter Brune if (flg) { 31107144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 31207144faaSPeter Brune } 313ee78dd50SPeter Brune 314fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 315ab8d36c9SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 316ab8d36c9SPeter Brune if (flg) { 317ab8d36c9SPeter Brune ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 318fde0ff24SPeter Brune } 319fde0ff24SPeter Brune 320ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 321ab8d36c9SPeter Brune if (flg) { 322ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 323ab8d36c9SPeter Brune } 324ee78dd50SPeter Brune 32566585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 326162d76ddSPeter Brune 32766585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 328162d76ddSPeter Brune 329c8c899caSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 330c8c899caSPeter Brune if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 331ab8d36c9SPeter Brune } 332ee78dd50SPeter Brune 333421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3348cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 335162d76ddSPeter Brune if (upflg) { 33666585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 337162d76ddSPeter Brune } 338162d76ddSPeter Brune if (downflg) { 33966585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 340162d76ddSPeter Brune } 341eff52c0eSPeter Brune 3429e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3439e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3449e764e56SPeter Brune if (!snes->linesearch) { 345f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 3461a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 3479e764e56SPeter Brune } 3489e764e56SPeter Brune } 3499e764e56SPeter Brune 350ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 351ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 352ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 353421d9b32SPeter Brune PetscFunctionReturn(0); 354421d9b32SPeter Brune } 355421d9b32SPeter Brune 356421d9b32SPeter Brune #undef __FUNCT__ 357421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 358421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 359421d9b32SPeter Brune { 360421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 361*656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 362ab8d36c9SPeter Brune PetscInt i; 363421d9b32SPeter Brune PetscErrorCode ierr; 364ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 365421d9b32SPeter Brune 366421d9b32SPeter Brune PetscFunctionBegin; 367ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 368ab8d36c9SPeter Brune if (isFine) { 369251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 370*656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 371421d9b32SPeter Brune if (iascii) { 372ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 373ab8d36c9SPeter Brune if (fas->galerkin) { 374ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 375421d9b32SPeter Brune } else { 376ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 377421d9b32SPeter Brune } 378ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 379ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 380ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 381ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 382ab8d36c9SPeter Brune if (!i) { 383ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 384421d9b32SPeter Brune } else { 385ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 386421d9b32SPeter Brune } 387ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 388ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 389ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 390ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 391ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 392ab8d36c9SPeter Brune } else if (i) { 393ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 394ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 395ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 396ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 397ab8d36c9SPeter Brune } 398ab8d36c9SPeter Brune } 399*656ede7eSPeter Brune } else if (isdraw) { 400*656ede7eSPeter Brune ierr = PetscPrintf(PETSC_COMM_WORLD,"trying to draw");CHKERRQ(ierr); 401*656ede7eSPeter Brune PetscDraw draw; 402*656ede7eSPeter Brune PetscReal x,y,bottom,th,wth; 403*656ede7eSPeter Brune SNES_FAS *curfas = fas; 404*656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 405*656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 406*656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 407*656ede7eSPeter Brune bottom = y - th; 408*656ede7eSPeter Brune while (curfas) { 409*656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 410*656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 411*656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 412*656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 413*656ede7eSPeter Brune bottom -= 5*th; 414*656ede7eSPeter Brune if (curfas->next) { 415*656ede7eSPeter Brune curfas = (SNES_FAS*)curfas->next->data; 416*656ede7eSPeter Brune } else { 417*656ede7eSPeter Brune curfas = PETSC_NULL; 418*656ede7eSPeter Brune } 419*656ede7eSPeter Brune } 420421d9b32SPeter Brune } else { 421421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 422421d9b32SPeter Brune } 423ab8d36c9SPeter Brune } 424421d9b32SPeter Brune PetscFunctionReturn(0); 425421d9b32SPeter Brune } 426421d9b32SPeter Brune 427421d9b32SPeter Brune #undef __FUNCT__ 42891f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private" 42939bd7f45SPeter Brune /* 43039bd7f45SPeter Brune Defines the action of the downsmoother 43139bd7f45SPeter Brune */ 43291f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 433b9c2fdf1SPeter Brune { 43439bd7f45SPeter Brune PetscErrorCode ierr = 0; 435742fe5e2SPeter Brune SNESConvergedReason reason; 436ab8d36c9SPeter Brune Vec FPC; 437ab8d36c9SPeter Brune SNES smoothd; 438421d9b32SPeter Brune PetscFunctionBegin; 439ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 440e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 441e4ed7901SPeter Brune ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr); 442ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 443742fe5e2SPeter Brune /* check convergence reason for the smoother */ 444ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 445e70c42e5SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) { 446742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 447742fe5e2SPeter Brune PetscFunctionReturn(0); 448742fe5e2SPeter Brune } 449ab8d36c9SPeter Brune ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 4504b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 451b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr); 45239bd7f45SPeter Brune PetscFunctionReturn(0); 45339bd7f45SPeter Brune } 45439bd7f45SPeter Brune 45539bd7f45SPeter Brune 45639bd7f45SPeter Brune #undef __FUNCT__ 45791f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private" 45839bd7f45SPeter Brune /* 45907144faaSPeter Brune Defines the action of the upsmoother 46039bd7f45SPeter Brune */ 46191f99d7cSPeter Brune PetscErrorCode SNESFASUpSmooth_Private (SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) { 46239bd7f45SPeter Brune PetscErrorCode ierr = 0; 46339bd7f45SPeter Brune SNESConvergedReason reason; 464ab8d36c9SPeter Brune Vec FPC; 465ab8d36c9SPeter Brune SNES smoothu; 46639bd7f45SPeter Brune PetscFunctionBegin; 467ab8d36c9SPeter Brune 468ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 469ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 47039bd7f45SPeter Brune /* check convergence reason for the smoother */ 471ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 47239bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 47339bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 47439bd7f45SPeter Brune PetscFunctionReturn(0); 47539bd7f45SPeter Brune } 476ab8d36c9SPeter Brune ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 4774b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 478b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr); 47939bd7f45SPeter Brune PetscFunctionReturn(0); 48039bd7f45SPeter Brune } 48139bd7f45SPeter Brune 48239bd7f45SPeter Brune #undef __FUNCT__ 483938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 484938e4a01SJed Brown /*@ 485938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 486938e4a01SJed Brown 487938e4a01SJed Brown Collective 488938e4a01SJed Brown 489938e4a01SJed Brown Input Arguments: 490938e4a01SJed Brown . snes - SNESFAS 491938e4a01SJed Brown 492938e4a01SJed Brown Output Arguments: 493938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 494938e4a01SJed Brown 495938e4a01SJed Brown Level: developer 496938e4a01SJed Brown 497938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 498938e4a01SJed Brown @*/ 499938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 500938e4a01SJed Brown { 501938e4a01SJed Brown PetscErrorCode ierr; 502938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 503938e4a01SJed Brown 504938e4a01SJed Brown PetscFunctionBegin; 505938e4a01SJed Brown if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 506938e4a01SJed Brown else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 507938e4a01SJed Brown else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 508938e4a01SJed Brown PetscFunctionReturn(0); 509938e4a01SJed Brown } 510938e4a01SJed Brown 511e9923e8dSJed Brown #undef __FUNCT__ 512e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 513e9923e8dSJed Brown /*@ 514e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 515e9923e8dSJed Brown 516e9923e8dSJed Brown Collective 517e9923e8dSJed Brown 518e9923e8dSJed Brown Input Arguments: 519e9923e8dSJed Brown + fine - SNES from which to restrict 520e9923e8dSJed Brown - Xfine - vector to restrict 521e9923e8dSJed Brown 522e9923e8dSJed Brown Output Arguments: 523e9923e8dSJed Brown . Xcoarse - result of restriction 524e9923e8dSJed Brown 525e9923e8dSJed Brown Level: developer 526e9923e8dSJed Brown 527e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 528e9923e8dSJed Brown @*/ 529e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 530e9923e8dSJed Brown { 531e9923e8dSJed Brown PetscErrorCode ierr; 532e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 533e9923e8dSJed Brown 534e9923e8dSJed Brown PetscFunctionBegin; 535e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 536e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 537e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 538e9923e8dSJed Brown if (fas->inject) { 539e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 540e9923e8dSJed Brown } else { 541e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 542e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 543e9923e8dSJed Brown } 544e9923e8dSJed Brown PetscFunctionReturn(0); 545e9923e8dSJed Brown } 546e9923e8dSJed Brown 547e9923e8dSJed Brown #undef __FUNCT__ 5488c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection" 54939bd7f45SPeter Brune /* 55039bd7f45SPeter Brune 55139bd7f45SPeter Brune Performs the FAS coarse correction as: 55239bd7f45SPeter Brune 55339bd7f45SPeter Brune fine problem: F(x) = 0 55439bd7f45SPeter Brune coarse problem: F^c(x) = b^c 55539bd7f45SPeter Brune 55639bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 55739bd7f45SPeter Brune 55839bd7f45SPeter Brune */ 5598c40d5fbSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 56039bd7f45SPeter Brune PetscErrorCode ierr; 56139bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 56239bd7f45SPeter Brune SNESConvergedReason reason; 563ab8d36c9SPeter Brune SNES next; 564ab8d36c9SPeter Brune Mat restrct, interpolate; 56539bd7f45SPeter Brune PetscFunctionBegin; 566ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 567ab8d36c9SPeter Brune if (next) { 568ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 569ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 570ab8d36c9SPeter Brune 571ab8d36c9SPeter Brune X_c = next->vec_sol; 572ab8d36c9SPeter Brune Xo_c = next->work[0]; 573ab8d36c9SPeter Brune F_c = next->vec_func; 574ab8d36c9SPeter Brune B_c = next->vec_rhs; 575efe1f98aSPeter Brune 576938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 577293a7e31SPeter Brune /* restrict the defect */ 578ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 579b9c2fdf1SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 580ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 581e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 582b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 583e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 584ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 585ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 586c90fad12SPeter Brune 587c90fad12SPeter Brune /* recurse to the next level */ 588e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 589ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 590ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 591742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 592742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 593742fe5e2SPeter Brune PetscFunctionReturn(0); 594742fe5e2SPeter Brune } 595fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 596fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 597ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 598293a7e31SPeter Brune } 59939bd7f45SPeter Brune PetscFunctionReturn(0); 60039bd7f45SPeter Brune } 60139bd7f45SPeter Brune 60239bd7f45SPeter Brune #undef __FUNCT__ 6032cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive" 60439bd7f45SPeter Brune /* 60539bd7f45SPeter Brune 60639bd7f45SPeter Brune The additive cycle looks like: 60739bd7f45SPeter Brune 60807144faaSPeter Brune xhat = x 60907144faaSPeter Brune xhat = dS(x, b) 61007144faaSPeter Brune x = coarsecorrection(xhat, b_d) 61107144faaSPeter Brune x = x + nu*(xhat - x); 61239bd7f45SPeter Brune (optional) x = uS(x, b) 61339bd7f45SPeter Brune 61439bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 61539bd7f45SPeter Brune 61639bd7f45SPeter Brune */ 61791f99d7cSPeter Brune PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) { 61807144faaSPeter Brune Vec F, B, Xhat; 61922c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 62039bd7f45SPeter Brune PetscErrorCode ierr; 62107144faaSPeter Brune SNESConvergedReason reason; 62222c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 62322c1e704SPeter Brune PetscBool lssuccess; 624ab8d36c9SPeter Brune SNES next; 625ab8d36c9SPeter Brune Mat restrct, interpolate; 62639bd7f45SPeter Brune PetscFunctionBegin; 627ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 62839bd7f45SPeter Brune F = snes->vec_func; 62939bd7f45SPeter Brune B = snes->vec_rhs; 630e7f468e7SPeter Brune Xhat = snes->work[1]; 63107144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 63207144faaSPeter Brune /* recurse first */ 633ab8d36c9SPeter Brune if (next) { 634ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 635ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 63607144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 637c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 638ab8d36c9SPeter Brune X_c = next->vec_sol; 639ab8d36c9SPeter Brune Xo_c = next->work[0]; 640ab8d36c9SPeter Brune F_c = next->vec_func; 641ab8d36c9SPeter Brune B_c = next->vec_rhs; 64239bd7f45SPeter Brune 643938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 64407144faaSPeter Brune /* restrict the defect */ 645ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 64607144faaSPeter Brune 64707144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 648ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 649e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 650b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 651e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 65207144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 65307144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 65407144faaSPeter Brune 65507144faaSPeter Brune /* recurse */ 656e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 657ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 65807144faaSPeter Brune 65907144faaSPeter Brune /* smooth on this level */ 66091f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 66107144faaSPeter Brune 662ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 66307144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 66407144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 66507144faaSPeter Brune PetscFunctionReturn(0); 66607144faaSPeter Brune } 66707144faaSPeter Brune 66807144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 669c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 670ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 67107144faaSPeter Brune 672ddebd997SPeter Brune /* additive correction of the coarse direction*/ 673f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 674f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 6759e764e56SPeter Brune if (!lssuccess) { 6769e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 6779e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 6789e764e56SPeter Brune PetscFunctionReturn(0); 6799e764e56SPeter Brune } 6809e764e56SPeter Brune } 681b9c2fdf1SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 68207144faaSPeter Brune } else { 68391f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 68407144faaSPeter Brune } 68539bd7f45SPeter Brune PetscFunctionReturn(0); 68639bd7f45SPeter Brune } 68739bd7f45SPeter Brune 68839bd7f45SPeter Brune #undef __FUNCT__ 6892cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative" 69039bd7f45SPeter Brune /* 69139bd7f45SPeter Brune 69239bd7f45SPeter Brune Defines the FAS cycle as: 69339bd7f45SPeter Brune 69439bd7f45SPeter Brune fine problem: F(x) = 0 69539bd7f45SPeter Brune coarse problem: F^c(x) = b^c 69639bd7f45SPeter Brune 69739bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 69839bd7f45SPeter Brune 69939bd7f45SPeter Brune correction: 70039bd7f45SPeter Brune 70139bd7f45SPeter Brune x = x + I(x^c - Rx) 70239bd7f45SPeter Brune 70339bd7f45SPeter Brune */ 70491f99d7cSPeter Brune PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) { 70539bd7f45SPeter Brune 70639bd7f45SPeter Brune PetscErrorCode ierr; 70739bd7f45SPeter Brune Vec F,B; 70839bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 70939bd7f45SPeter Brune 71039bd7f45SPeter Brune PetscFunctionBegin; 71139bd7f45SPeter Brune F = snes->vec_func; 71239bd7f45SPeter Brune B = snes->vec_rhs; 71339bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 71491f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 715c90fad12SPeter Brune if (fas->level != 0) { 7168c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 71791f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 718fe6f9142SPeter Brune } 719fa9694d7SPeter Brune PetscFunctionReturn(0); 720421d9b32SPeter Brune } 721421d9b32SPeter Brune 722421d9b32SPeter Brune #undef __FUNCT__ 723421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 724421d9b32SPeter Brune 725421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 726421d9b32SPeter Brune { 727fa9694d7SPeter Brune PetscErrorCode ierr; 728fe6f9142SPeter Brune PetscInt i, maxits; 729ddb5aff1SPeter Brune Vec X, F; 730fe6f9142SPeter Brune PetscReal fnorm; 731b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 732b17ce1afSJed Brown DM dm; 733e70c42e5SPeter Brune PetscBool isFine; 734b17ce1afSJed Brown 735421d9b32SPeter Brune PetscFunctionBegin; 736fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 737fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 738fa9694d7SPeter Brune X = snes->vec_sol; 739f5a6d4f9SBarry Smith F = snes->vec_func; 740293a7e31SPeter Brune 74118a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 742293a7e31SPeter Brune /*norm setup */ 743fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 744fe6f9142SPeter Brune snes->iter = 0; 745fe6f9142SPeter Brune snes->norm = 0.; 746fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 747e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 748fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 749fe6f9142SPeter Brune if (snes->domainerror) { 750fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 751fe6f9142SPeter Brune PetscFunctionReturn(0); 752fe6f9142SPeter Brune } 753e4ed7901SPeter Brune } else { 754e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 755e4ed7901SPeter Brune } 756e4ed7901SPeter Brune 757e4ed7901SPeter Brune if (!snes->norm_init_set) { 758fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 759fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 760fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 761e4ed7901SPeter Brune } else { 762e4ed7901SPeter Brune fnorm = snes->norm_init; 763e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 764e4ed7901SPeter Brune } 765e4ed7901SPeter Brune 766fe6f9142SPeter Brune snes->norm = fnorm; 767fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 768fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 769fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 770fe6f9142SPeter Brune 771fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 772fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 773fe6f9142SPeter Brune /* test convergence */ 774fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 775fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 776e4ed7901SPeter Brune 777b17ce1afSJed Brown 778b9c2fdf1SPeter Brune if (isFine) { 779b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 780b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 781b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 782b17ce1afSJed Brown DM dmcoarse; 783b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 784b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 785b17ce1afSJed Brown dm = dmcoarse; 786b17ce1afSJed Brown } 787b9c2fdf1SPeter Brune } 788b17ce1afSJed Brown 789fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 790fe6f9142SPeter Brune /* Call general purpose update function */ 791646217ecSPeter Brune 792fe6f9142SPeter Brune if (snes->ops->update) { 793fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 794fe6f9142SPeter Brune } 79507144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 79691f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 79707144faaSPeter Brune } else { 79891f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 79907144faaSPeter Brune } 800742fe5e2SPeter Brune 801742fe5e2SPeter Brune /* check for FAS cycle divergence */ 802742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 803742fe5e2SPeter Brune PetscFunctionReturn(0); 804742fe5e2SPeter Brune } 805b9c2fdf1SPeter Brune 806c90fad12SPeter Brune /* Monitor convergence */ 807c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 808c90fad12SPeter Brune snes->iter = i+1; 809c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 810c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 811c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 812c90fad12SPeter Brune /* Test for convergence */ 81366585501SPeter Brune if (isFine) { 814b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 815c90fad12SPeter Brune if (snes->reason) break; 816fe6f9142SPeter Brune } 81766585501SPeter Brune } 818fe6f9142SPeter Brune if (i == maxits) { 819fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 820fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 821fe6f9142SPeter Brune } 822421d9b32SPeter Brune PetscFunctionReturn(0); 823421d9b32SPeter Brune } 824