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; 1406b2b7091SBarry Smith PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 1416b2b7091SBarry Smith PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool *,PetscBool *,void*); 142eff52c0eSPeter Brune 1436b2b7091SBarry Smith PetscFunctionBegin; 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 } 182ab8d36c9SPeter Brune 18379d9a41aSPeter Brune if (snes->dm) { 184ab8d36c9SPeter Brune /* set the smoother DMs properly */ 185ab8d36c9SPeter Brune if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 186ab8d36c9SPeter Brune ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 18779d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 188ab8d36c9SPeter Brune if (next) { 18979d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 190ab8d36c9SPeter Brune if (!next->dm) { 191ab8d36c9SPeter Brune ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr); 192ab8d36c9SPeter Brune ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 19379d9a41aSPeter Brune } 19479d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 19579d9a41aSPeter Brune if (!fas->interpolate) { 196ab8d36c9SPeter Brune ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 197bccf9bb3SJed Brown if (!fas->restrct) { 198bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 19979d9a41aSPeter Brune fas->restrct = fas->interpolate; 20079d9a41aSPeter Brune } 201bccf9bb3SJed Brown } 20279d9a41aSPeter Brune /* set the injection from the DM */ 20379d9a41aSPeter Brune if (!fas->inject) { 204ab8d36c9SPeter Brune ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 20579d9a41aSPeter Brune ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr); 20679d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 20779d9a41aSPeter Brune } 20879d9a41aSPeter Brune } 20979d9a41aSPeter Brune } 21079d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 21179d9a41aSPeter Brune if (fas->galerkin) { 212ab8d36c9SPeter Brune if (next) 213ab8d36c9SPeter Brune ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 214ab8d36c9SPeter Brune if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 215ab8d36c9SPeter Brune if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 21679d9a41aSPeter Brune } 21779d9a41aSPeter Brune 218534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 219534ebe21SPeter Brune if (fas->smoothd) { 220bc3f2f05SPeter Brune if (fas->level == 0 && fas->levels != 1) { 221534ebe21SPeter Brune ierr = SNESSetNormType(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); 222534ebe21SPeter Brune } else { 223534ebe21SPeter Brune ierr = SNESSetNormType(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 224534ebe21SPeter Brune } 2257fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr); 226534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); 227f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 228f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr); 2296b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2306b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2316b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2326b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 233f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 234534ebe21SPeter Brune } 235534ebe21SPeter Brune 236534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 237534ebe21SPeter Brune if (fas->smoothu) { 238534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 239534ebe21SPeter Brune ierr = SNESSetNormType(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); 240534ebe21SPeter Brune } else { 241534ebe21SPeter Brune ierr = SNESSetNormType(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 242534ebe21SPeter Brune } 2437fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr); 244534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); 245f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 246f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr); 2476b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2486b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2496b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2506b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 251f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 252534ebe21SPeter Brune } 253d06165b7SPeter Brune 254ab8d36c9SPeter Brune if (next) { 25579d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 256ab8d36c9SPeter Brune if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 257ab8d36c9SPeter Brune if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 2587fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); 259f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr); 260f89ba88eSPeter Brune ierr = SNESGetSNESLineSearch(fas->next,&slinesearch);CHKERRQ(ierr); 2616b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2626b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2636b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2646b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 265f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 266ab8d36c9SPeter Brune ierr = SNESSetUp(next);CHKERRQ(ierr); 26779d9a41aSPeter Brune } 2686273346dSPeter Brune /* setup FAS work vectors */ 2696273346dSPeter Brune if (fas->galerkin) { 2706273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 2716273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 2726273346dSPeter Brune } 273421d9b32SPeter Brune PetscFunctionReturn(0); 274421d9b32SPeter Brune } 275421d9b32SPeter Brune 276421d9b32SPeter Brune #undef __FUNCT__ 277421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 278421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 279421d9b32SPeter Brune { 280ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 281ee78dd50SPeter Brune PetscInt levels = 1; 2824d26bfa5SPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE; 283421d9b32SPeter Brune PetscErrorCode ierr; 284ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 28507144faaSPeter Brune SNESFASType fastype; 286fde0ff24SPeter Brune const char *optionsprefix; 287f1c6b773SPeter Brune SNESLineSearch linesearch; 28866585501SPeter Brune PetscInt m, n_up, n_down; 289ab8d36c9SPeter Brune SNES next; 290ab8d36c9SPeter Brune PetscBool isFine; 291421d9b32SPeter Brune 292421d9b32SPeter Brune PetscFunctionBegin; 293ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 294c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 295ee78dd50SPeter Brune 296ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 297ab8d36c9SPeter Brune if (isFine) { 298ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 299c732cbdbSBarry Smith if (!flg && snes->dm) { 300c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 301c732cbdbSBarry Smith levels++; 302d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 303c732cbdbSBarry Smith } 304ee78dd50SPeter Brune ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr); 30507144faaSPeter Brune fastype = fas->fastype; 30607144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 30707144faaSPeter Brune if (flg) { 30807144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 30907144faaSPeter Brune } 310ee78dd50SPeter Brune 311fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 312ab8d36c9SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 313ab8d36c9SPeter Brune if (flg) { 314ab8d36c9SPeter Brune ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 315fde0ff24SPeter Brune } 316fde0ff24SPeter Brune 317ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 318ab8d36c9SPeter Brune if (flg) { 319ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 320ab8d36c9SPeter Brune } 321ee78dd50SPeter Brune 32266585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 323162d76ddSPeter Brune 32466585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 325162d76ddSPeter Brune 326c8c899caSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 327c8c899caSPeter Brune if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 328ab8d36c9SPeter Brune } 329ee78dd50SPeter Brune 330421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3318cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 332162d76ddSPeter Brune if (upflg) { 33366585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 334162d76ddSPeter Brune } 335162d76ddSPeter Brune if (downflg) { 33666585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 337162d76ddSPeter Brune } 338eff52c0eSPeter Brune 3399e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3409e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3419e764e56SPeter Brune if (!snes->linesearch) { 342f1c6b773SPeter Brune ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr); 3431a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 3449e764e56SPeter Brune } 3459e764e56SPeter Brune } 3469e764e56SPeter Brune 347ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 348ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 349ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 350421d9b32SPeter Brune PetscFunctionReturn(0); 351421d9b32SPeter Brune } 352421d9b32SPeter Brune 353421d9b32SPeter Brune #undef __FUNCT__ 354421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 355421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 356421d9b32SPeter Brune { 357421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *) snes->data; 358656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 359ab8d36c9SPeter Brune PetscInt i; 360421d9b32SPeter Brune PetscErrorCode ierr; 361ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 362421d9b32SPeter Brune 363421d9b32SPeter Brune PetscFunctionBegin; 364ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 365ab8d36c9SPeter Brune if (isFine) { 366251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 367656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 368421d9b32SPeter Brune if (iascii) { 369ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 370ab8d36c9SPeter Brune if (fas->galerkin) { 371ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 372421d9b32SPeter Brune } else { 373ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 374421d9b32SPeter Brune } 375ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 376ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 377ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 378ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 379ab8d36c9SPeter Brune if (!i) { 380ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 381421d9b32SPeter Brune } else { 382ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 383421d9b32SPeter Brune } 384ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 385ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 386ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 387ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 388ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 389ab8d36c9SPeter Brune } else if (i) { 390ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 391ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 392ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 393ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 394ab8d36c9SPeter Brune } 395ab8d36c9SPeter Brune } 396656ede7eSPeter Brune } else if (isdraw) { 397656ede7eSPeter Brune PetscDraw draw; 398b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 399656ede7eSPeter Brune SNES_FAS *curfas = fas; 400656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 401656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 402656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 403656ede7eSPeter Brune bottom = y - th; 404656ede7eSPeter Brune while (curfas) { 405b4375e8dSPeter Brune if (!curfas->smoothu) { 406656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 407656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 408656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 409b4375e8dSPeter Brune } else { 410b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 411b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 412b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 413b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 414b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 415b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 416b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 417b4375e8dSPeter Brune } 418656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 419656ede7eSPeter Brune bottom -= 5*th; 420656ede7eSPeter Brune if (curfas->next) { 421656ede7eSPeter Brune curfas = (SNES_FAS*)curfas->next->data; 422656ede7eSPeter Brune } else { 423656ede7eSPeter Brune curfas = PETSC_NULL; 424656ede7eSPeter Brune } 425656ede7eSPeter Brune } 426421d9b32SPeter Brune } else { 427421d9b32SPeter Brune SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name); 428421d9b32SPeter Brune } 429ab8d36c9SPeter Brune } 430421d9b32SPeter Brune PetscFunctionReturn(0); 431421d9b32SPeter Brune } 432421d9b32SPeter Brune 433421d9b32SPeter Brune #undef __FUNCT__ 43491f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private" 43539bd7f45SPeter Brune /* 43639bd7f45SPeter Brune Defines the action of the downsmoother 43739bd7f45SPeter Brune */ 43891f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 439b9c2fdf1SPeter Brune { 44039bd7f45SPeter Brune PetscErrorCode ierr = 0; 441742fe5e2SPeter Brune SNESConvergedReason reason; 442ab8d36c9SPeter Brune Vec FPC; 443ab8d36c9SPeter Brune SNES smoothd; 444*6e111a19SKarl Rupp 445421d9b32SPeter Brune PetscFunctionBegin; 446ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 447e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 448e4ed7901SPeter Brune ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr); 449ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 450742fe5e2SPeter Brune /* check convergence reason for the smoother */ 451ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 452e70c42e5SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) { 453742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 454742fe5e2SPeter Brune PetscFunctionReturn(0); 455742fe5e2SPeter Brune } 456ab8d36c9SPeter Brune ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 4574b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 458b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr); 45939bd7f45SPeter Brune PetscFunctionReturn(0); 46039bd7f45SPeter Brune } 46139bd7f45SPeter Brune 46239bd7f45SPeter Brune 46339bd7f45SPeter Brune #undef __FUNCT__ 46491f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private" 46539bd7f45SPeter Brune /* 46607144faaSPeter Brune Defines the action of the upsmoother 46739bd7f45SPeter Brune */ 4680adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private (SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 4690adebc6cSBarry Smith { 47039bd7f45SPeter Brune PetscErrorCode ierr = 0; 47139bd7f45SPeter Brune SNESConvergedReason reason; 472ab8d36c9SPeter Brune Vec FPC; 473ab8d36c9SPeter Brune SNES smoothu; 474ab8d36c9SPeter Brune 475*6e111a19SKarl Rupp PetscFunctionBegin; 476ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 477ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 47839bd7f45SPeter Brune /* check convergence reason for the smoother */ 479ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 48039bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 48139bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 48239bd7f45SPeter Brune PetscFunctionReturn(0); 48339bd7f45SPeter Brune } 484ab8d36c9SPeter Brune ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr); 4854b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 486b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr); 48739bd7f45SPeter Brune PetscFunctionReturn(0); 48839bd7f45SPeter Brune } 48939bd7f45SPeter Brune 49039bd7f45SPeter Brune #undef __FUNCT__ 491938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 492938e4a01SJed Brown /*@ 493938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 494938e4a01SJed Brown 495938e4a01SJed Brown Collective 496938e4a01SJed Brown 497938e4a01SJed Brown Input Arguments: 498938e4a01SJed Brown . snes - SNESFAS 499938e4a01SJed Brown 500938e4a01SJed Brown Output Arguments: 501938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 502938e4a01SJed Brown 503938e4a01SJed Brown Level: developer 504938e4a01SJed Brown 505938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 506938e4a01SJed Brown @*/ 507938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 508938e4a01SJed Brown { 509938e4a01SJed Brown PetscErrorCode ierr; 510938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 511938e4a01SJed Brown 512938e4a01SJed Brown PetscFunctionBegin; 513938e4a01SJed Brown if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);} 514938e4a01SJed Brown else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);} 515938e4a01SJed Brown else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 516938e4a01SJed Brown PetscFunctionReturn(0); 517938e4a01SJed Brown } 518938e4a01SJed Brown 519e9923e8dSJed Brown #undef __FUNCT__ 520e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 521e9923e8dSJed Brown /*@ 522e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 523e9923e8dSJed Brown 524e9923e8dSJed Brown Collective 525e9923e8dSJed Brown 526e9923e8dSJed Brown Input Arguments: 527e9923e8dSJed Brown + fine - SNES from which to restrict 528e9923e8dSJed Brown - Xfine - vector to restrict 529e9923e8dSJed Brown 530e9923e8dSJed Brown Output Arguments: 531e9923e8dSJed Brown . Xcoarse - result of restriction 532e9923e8dSJed Brown 533e9923e8dSJed Brown Level: developer 534e9923e8dSJed Brown 535e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 536e9923e8dSJed Brown @*/ 537e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 538e9923e8dSJed Brown { 539e9923e8dSJed Brown PetscErrorCode ierr; 540e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 541e9923e8dSJed Brown 542e9923e8dSJed Brown PetscFunctionBegin; 543e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 544e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 545e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 546e9923e8dSJed Brown if (fas->inject) { 547e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 548e9923e8dSJed Brown } else { 549e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 550e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 551e9923e8dSJed Brown } 552e9923e8dSJed Brown PetscFunctionReturn(0); 553e9923e8dSJed Brown } 554e9923e8dSJed Brown 555e9923e8dSJed Brown #undef __FUNCT__ 5568c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection" 55739bd7f45SPeter Brune /* 55839bd7f45SPeter Brune 55939bd7f45SPeter Brune Performs the FAS coarse correction as: 56039bd7f45SPeter Brune 56139bd7f45SPeter Brune fine problem: F(x) = 0 56239bd7f45SPeter Brune coarse problem: F^c(x) = b^c 56339bd7f45SPeter Brune 56439bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 56539bd7f45SPeter Brune 56639bd7f45SPeter Brune */ 5670adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 5680adebc6cSBarry Smith { 56939bd7f45SPeter Brune PetscErrorCode ierr; 57039bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 57139bd7f45SPeter Brune SNESConvergedReason reason; 572ab8d36c9SPeter Brune SNES next; 573ab8d36c9SPeter Brune Mat restrct, interpolate; 57439bd7f45SPeter Brune PetscFunctionBegin; 575ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 576ab8d36c9SPeter Brune if (next) { 577ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 578ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 579ab8d36c9SPeter Brune 580ab8d36c9SPeter Brune X_c = next->vec_sol; 581ab8d36c9SPeter Brune Xo_c = next->work[0]; 582ab8d36c9SPeter Brune F_c = next->vec_func; 583ab8d36c9SPeter Brune B_c = next->vec_rhs; 584efe1f98aSPeter Brune 585938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 586293a7e31SPeter Brune /* restrict the defect */ 587ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 588b9c2fdf1SPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 589ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 590e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 591b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 592e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 593ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 594ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 595c90fad12SPeter Brune 596c90fad12SPeter Brune /* recurse to the next level */ 597e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 598ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 599ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 600742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 601742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 602742fe5e2SPeter Brune PetscFunctionReturn(0); 603742fe5e2SPeter Brune } 604fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 605fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 606ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 607293a7e31SPeter Brune } 60839bd7f45SPeter Brune PetscFunctionReturn(0); 60939bd7f45SPeter Brune } 61039bd7f45SPeter Brune 61139bd7f45SPeter Brune #undef __FUNCT__ 6122cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive" 61339bd7f45SPeter Brune /* 61439bd7f45SPeter Brune 61539bd7f45SPeter Brune The additive cycle looks like: 61639bd7f45SPeter Brune 61707144faaSPeter Brune xhat = x 61807144faaSPeter Brune xhat = dS(x, b) 61907144faaSPeter Brune x = coarsecorrection(xhat, b_d) 62007144faaSPeter Brune x = x + nu*(xhat - x); 62139bd7f45SPeter Brune (optional) x = uS(x, b) 62239bd7f45SPeter Brune 62339bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 62439bd7f45SPeter Brune 62539bd7f45SPeter Brune */ 6260adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 6270adebc6cSBarry Smith { 62807144faaSPeter Brune Vec F, B, Xhat; 62922c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 63039bd7f45SPeter Brune PetscErrorCode ierr; 63107144faaSPeter Brune SNESConvergedReason reason; 63222c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 63322c1e704SPeter Brune PetscBool lssuccess; 634ab8d36c9SPeter Brune SNES next; 635ab8d36c9SPeter Brune Mat restrct, interpolate; 63639bd7f45SPeter Brune PetscFunctionBegin; 637ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 63839bd7f45SPeter Brune F = snes->vec_func; 63939bd7f45SPeter Brune B = snes->vec_rhs; 640e7f468e7SPeter Brune Xhat = snes->work[1]; 64107144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 64207144faaSPeter Brune /* recurse first */ 643ab8d36c9SPeter Brune if (next) { 644ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 645ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 64607144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 647c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 648ab8d36c9SPeter Brune X_c = next->vec_sol; 649ab8d36c9SPeter Brune Xo_c = next->work[0]; 650ab8d36c9SPeter Brune F_c = next->vec_func; 651ab8d36c9SPeter Brune B_c = next->vec_rhs; 65239bd7f45SPeter Brune 653938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 65407144faaSPeter Brune /* restrict the defect */ 655ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 65607144faaSPeter Brune 65707144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 658ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 659e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 660b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 661e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 66207144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 66307144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 66407144faaSPeter Brune 66507144faaSPeter Brune /* recurse */ 666e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 667ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 66807144faaSPeter Brune 66907144faaSPeter Brune /* smooth on this level */ 67091f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 67107144faaSPeter Brune 672ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 67307144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 67407144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 67507144faaSPeter Brune PetscFunctionReturn(0); 67607144faaSPeter Brune } 67707144faaSPeter Brune 67807144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 679c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 680ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 68107144faaSPeter Brune 682ddebd997SPeter Brune /* additive correction of the coarse direction*/ 683f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 684f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 6859e764e56SPeter Brune if (!lssuccess) { 6869e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 6879e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 6889e764e56SPeter Brune PetscFunctionReturn(0); 6899e764e56SPeter Brune } 6909e764e56SPeter Brune } 691b9c2fdf1SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 69207144faaSPeter Brune } else { 69391f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 69407144faaSPeter Brune } 69539bd7f45SPeter Brune PetscFunctionReturn(0); 69639bd7f45SPeter Brune } 69739bd7f45SPeter Brune 69839bd7f45SPeter Brune #undef __FUNCT__ 6992cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative" 70039bd7f45SPeter Brune /* 70139bd7f45SPeter Brune 70239bd7f45SPeter Brune Defines the FAS cycle as: 70339bd7f45SPeter Brune 70439bd7f45SPeter Brune fine problem: F(x) = 0 70539bd7f45SPeter Brune coarse problem: F^c(x) = b^c 70639bd7f45SPeter Brune 70739bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 70839bd7f45SPeter Brune 70939bd7f45SPeter Brune correction: 71039bd7f45SPeter Brune 71139bd7f45SPeter Brune x = x + I(x^c - Rx) 71239bd7f45SPeter Brune 71339bd7f45SPeter Brune */ 7140adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7150adebc6cSBarry Smith { 71639bd7f45SPeter Brune 71739bd7f45SPeter Brune PetscErrorCode ierr; 71839bd7f45SPeter Brune Vec F,B; 71939bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 72039bd7f45SPeter Brune 72139bd7f45SPeter Brune PetscFunctionBegin; 72239bd7f45SPeter Brune F = snes->vec_func; 72339bd7f45SPeter Brune B = snes->vec_rhs; 72439bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 72591f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 726c90fad12SPeter Brune if (fas->level != 0) { 7278c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 72891f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 729fe6f9142SPeter Brune } 730fa9694d7SPeter Brune PetscFunctionReturn(0); 731421d9b32SPeter Brune } 732421d9b32SPeter Brune 733421d9b32SPeter Brune #undef __FUNCT__ 734421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 735421d9b32SPeter Brune 736421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 737421d9b32SPeter Brune { 738fa9694d7SPeter Brune PetscErrorCode ierr; 739fe6f9142SPeter Brune PetscInt i, maxits; 740ddb5aff1SPeter Brune Vec X, F; 741fe6f9142SPeter Brune PetscReal fnorm; 742b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS *)snes->data,*ffas; 743b17ce1afSJed Brown DM dm; 744e70c42e5SPeter Brune PetscBool isFine; 745b17ce1afSJed Brown 746421d9b32SPeter Brune PetscFunctionBegin; 747fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 748fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 749fa9694d7SPeter Brune X = snes->vec_sol; 750f5a6d4f9SBarry Smith F = snes->vec_func; 751293a7e31SPeter Brune 75218a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 753293a7e31SPeter Brune /*norm setup */ 754fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 755fe6f9142SPeter Brune snes->iter = 0; 756fe6f9142SPeter Brune snes->norm = 0.; 757fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 758e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 759fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 760fe6f9142SPeter Brune if (snes->domainerror) { 761fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 762fe6f9142SPeter Brune PetscFunctionReturn(0); 763fe6f9142SPeter Brune } 764e4ed7901SPeter Brune } else { 765e4ed7901SPeter Brune snes->vec_func_init_set = PETSC_FALSE; 766e4ed7901SPeter Brune } 767e4ed7901SPeter Brune 768e4ed7901SPeter Brune if (!snes->norm_init_set) { 769fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 770fe6f9142SPeter Brune if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm"); 771fe6f9142SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 772e4ed7901SPeter Brune } else { 773e4ed7901SPeter Brune fnorm = snes->norm_init; 774e4ed7901SPeter Brune snes->norm_init_set = PETSC_FALSE; 775e4ed7901SPeter Brune } 776e4ed7901SPeter Brune 777fe6f9142SPeter Brune snes->norm = fnorm; 778fe6f9142SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 779fe6f9142SPeter Brune SNESLogConvHistory(snes,fnorm,0); 780fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 781fe6f9142SPeter Brune 782fe6f9142SPeter Brune /* set parameter for default relative tolerance convergence test */ 783fe6f9142SPeter Brune snes->ttol = fnorm*snes->rtol; 784fe6f9142SPeter Brune /* test convergence */ 785fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 786fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 787e4ed7901SPeter Brune 788b17ce1afSJed Brown 789b9c2fdf1SPeter Brune if (isFine) { 790b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 791b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 792b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 793b17ce1afSJed Brown DM dmcoarse; 794b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 795b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 796b17ce1afSJed Brown dm = dmcoarse; 797b17ce1afSJed Brown } 798b9c2fdf1SPeter Brune } 799b17ce1afSJed Brown 800fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 801fe6f9142SPeter Brune /* Call general purpose update function */ 802646217ecSPeter Brune 803fe6f9142SPeter Brune if (snes->ops->update) { 804fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 805fe6f9142SPeter Brune } 80607144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 80791f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 80807144faaSPeter Brune } else { 80991f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 81007144faaSPeter Brune } 811742fe5e2SPeter Brune 812742fe5e2SPeter Brune /* check for FAS cycle divergence */ 813742fe5e2SPeter Brune if (snes->reason != SNES_CONVERGED_ITERATING) { 814742fe5e2SPeter Brune PetscFunctionReturn(0); 815742fe5e2SPeter Brune } 816b9c2fdf1SPeter Brune 817c90fad12SPeter Brune /* Monitor convergence */ 818c90fad12SPeter Brune ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr); 819c90fad12SPeter Brune snes->iter = i+1; 820c90fad12SPeter Brune ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr); 821c90fad12SPeter Brune SNESLogConvHistory(snes,snes->norm,0); 822c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 823c90fad12SPeter Brune /* Test for convergence */ 82466585501SPeter Brune if (isFine) { 825b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 826c90fad12SPeter Brune if (snes->reason) break; 827fe6f9142SPeter Brune } 82866585501SPeter Brune } 829fe6f9142SPeter Brune if (i == maxits) { 830fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 831fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 832fe6f9142SPeter Brune } 833421d9b32SPeter Brune PetscFunctionReturn(0); 834421d9b32SPeter Brune } 835