1421d9b32SPeter Brune /* Defines the basic SNES object */ 26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3421d9b32SPeter Brune 434d65b3cSPeter Brune const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","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 26928e959bSPeter Brune . -snes_fas_type<additive,multiplicative,full,kaskade> - 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 31928e959bSPeter Brune . -snes_fas_full_downsweepsmooth<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS 32d3bc2379SPeter Brune . -fas_levels_snes_ - SNES options for all smoothers 337d84e935SPeter Brune . -fas_levels_cycle_snes_ - SNES options for all cycles 34d3bc2379SPeter Brune . -fas_levels_i_snes_ - SNES options for the smoothers on level i 357d84e935SPeter Brune . -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i 36d3bc2379SPeter Brune - -fas_coarse_snes_ - SNES options for the coarsest smoother 37d3bc2379SPeter Brune 38d3bc2379SPeter Brune Notes: 39d3bc2379SPeter Brune The organization of the FAS solver is slightly different from the organization of PCMG 40d3bc2379SPeter Brune As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance. 41d3bc2379SPeter Brune The cycle SNES instance may be used for monitoring convergence on a particular level. 421fbfccc6SJed Brown 437d84e935SPeter Brune Level: beginner 441fbfccc6SJed Brown 45d3bc2379SPeter Brune .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 461fbfccc6SJed Brown M*/ 47421d9b32SPeter Brune 48421d9b32SPeter Brune #undef __FUNCT__ 49421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 508cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes) 51421d9b32SPeter Brune { 52421d9b32SPeter Brune SNES_FAS *fas; 53421d9b32SPeter Brune PetscErrorCode ierr; 54421d9b32SPeter Brune 55421d9b32SPeter Brune PetscFunctionBegin; 56421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 57421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 58421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 59421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 60421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 61421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 62421d9b32SPeter Brune 63ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 64ed020824SBarry Smith snes->usespc = PETSC_FALSE; 65ed020824SBarry Smith 6688976e71SPeter Brune if (!snes->tolerancesset) { 670e444f03SPeter Brune snes->max_funcs = 30000; 680e444f03SPeter Brune snes->max_its = 10000; 6988976e71SPeter Brune } 700e444f03SPeter Brune 71421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 721aa26658SKarl Rupp 73421d9b32SPeter Brune snes->data = (void*) fas; 74421d9b32SPeter Brune fas->level = 0; 75293a7e31SPeter Brune fas->levels = 1; 76ee78dd50SPeter Brune fas->n_cycles = 1; 77ee78dd50SPeter Brune fas->max_up_it = 1; 78ee78dd50SPeter Brune fas->max_down_it = 1; 790298fd71SBarry Smith fas->smoothu = NULL; 800298fd71SBarry Smith fas->smoothd = NULL; 810298fd71SBarry Smith fas->next = NULL; 820298fd71SBarry Smith fas->previous = NULL; 83ab8d36c9SPeter Brune fas->fine = snes; 840298fd71SBarry Smith fas->interpolate = NULL; 850298fd71SBarry Smith fas->restrct = NULL; 860298fd71SBarry Smith fas->inject = NULL; 870298fd71SBarry Smith fas->monitor = NULL; 88cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 89ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 90928e959bSPeter Brune fas->full_downsweep = PETSC_FALSE; 910dd27c6cSPeter Brune 920298fd71SBarry Smith fas->eventsmoothsetup = 0; 930298fd71SBarry Smith fas->eventsmoothsolve = 0; 940298fd71SBarry Smith fas->eventresidual = 0; 950298fd71SBarry Smith fas->eventinterprestrict = 0; 96efe1f98aSPeter Brune PetscFunctionReturn(0); 97efe1f98aSPeter Brune } 98efe1f98aSPeter Brune 99421d9b32SPeter Brune #undef __FUNCT__ 100421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 101421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 102421d9b32SPeter Brune { 10377df8cc4SPeter Brune PetscErrorCode ierr = 0; 104421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 105421d9b32SPeter Brune 106421d9b32SPeter Brune PetscFunctionBegin; 107ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr); 108ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 1093dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 110bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 111bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 112bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 113742fe5e2SPeter Brune if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr); 114421d9b32SPeter Brune PetscFunctionReturn(0); 115421d9b32SPeter Brune } 116421d9b32SPeter Brune 117421d9b32SPeter Brune #undef __FUNCT__ 118421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 119421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 120421d9b32SPeter Brune { 121421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 122742fe5e2SPeter Brune PetscErrorCode ierr = 0; 123421d9b32SPeter Brune 124421d9b32SPeter Brune PetscFunctionBegin; 125421d9b32SPeter Brune /* recursively resets and then destroys */ 12679d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 1271aa26658SKarl Rupp if (fas->next) { 1281aa26658SKarl Rupp ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 1291aa26658SKarl Rupp } 130421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 131421d9b32SPeter Brune PetscFunctionReturn(0); 132421d9b32SPeter Brune } 133421d9b32SPeter Brune 134421d9b32SPeter Brune #undef __FUNCT__ 135421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 136421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 137421d9b32SPeter Brune { 13848bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 139421d9b32SPeter Brune PetscErrorCode ierr; 140efe1f98aSPeter Brune VecScatter injscatter; 141d1adcc6fSPeter Brune PetscInt dm_levels; 1423dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 143ab8d36c9SPeter Brune SNES next; 144ab8d36c9SPeter Brune PetscBool isFine; 145f89ba88eSPeter Brune SNESLineSearch linesearch; 146f89ba88eSPeter Brune SNESLineSearch slinesearch; 147f89ba88eSPeter Brune void *lsprectx,*lspostctx; 1486b2b7091SBarry Smith PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 1496b2b7091SBarry Smith PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 150eff52c0eSPeter Brune 1516b2b7091SBarry Smith PetscFunctionBegin; 152ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 153ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 154d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 155d1adcc6fSPeter Brune dm_levels++; 156cc05f883SPeter Brune if (dm_levels > fas->levels) { 1572e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 1583dccd265SPeter Brune vec_sol = snes->vec_sol; 1593dccd265SPeter Brune vec_func = snes->vec_func; 1603dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 1613dccd265SPeter Brune vec_rhs = snes->vec_rhs; 1620298fd71SBarry Smith snes->vec_sol = NULL; 1630298fd71SBarry Smith snes->vec_func = NULL; 1640298fd71SBarry Smith snes->vec_sol_update = NULL; 1650298fd71SBarry Smith snes->vec_rhs = NULL; 1663dccd265SPeter Brune 1673dccd265SPeter Brune /* reset the number of levels */ 1680298fd71SBarry Smith ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr); 169cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 1703dccd265SPeter Brune 1713dccd265SPeter Brune snes->vec_sol = vec_sol; 1723dccd265SPeter Brune snes->vec_func = vec_func; 1733dccd265SPeter Brune snes->vec_rhs = vec_rhs; 1743dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 175d1adcc6fSPeter Brune } 176d1adcc6fSPeter Brune } 177ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 178ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 1793dccd265SPeter Brune 180fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 181cc05f883SPeter Brune 182ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 183ab8d36c9SPeter Brune if (!fas->smoothd) { 184ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 185ab8d36c9SPeter Brune } 186ab8d36c9SPeter Brune 18779d9a41aSPeter Brune if (snes->dm) { 188ab8d36c9SPeter Brune /* set the smoother DMs properly */ 189ab8d36c9SPeter Brune if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 190ab8d36c9SPeter Brune ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 19179d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 192ab8d36c9SPeter Brune if (next) { 19379d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 194ab8d36c9SPeter Brune if (!next->dm) { 195ce94432eSBarry Smith ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr); 196ab8d36c9SPeter Brune ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 19779d9a41aSPeter Brune } 19879d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 19979d9a41aSPeter Brune if (!fas->interpolate) { 200ab8d36c9SPeter Brune ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 201bccf9bb3SJed Brown if (!fas->restrct) { 202bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 20379d9a41aSPeter Brune fas->restrct = fas->interpolate; 20479d9a41aSPeter Brune } 205bccf9bb3SJed Brown } 20679d9a41aSPeter Brune /* set the injection from the DM */ 20779d9a41aSPeter Brune if (!fas->inject) { 208ab8d36c9SPeter Brune ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 209ce94432eSBarry Smith ierr = MatCreateScatter(PetscObjectComm((PetscObject)snes), injscatter, &fas->inject);CHKERRQ(ierr); 21079d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr); 21179d9a41aSPeter Brune } 21279d9a41aSPeter Brune } 21379d9a41aSPeter Brune } 21479d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 21579d9a41aSPeter Brune if (fas->galerkin) { 2161aa26658SKarl Rupp if (next) { 2170298fd71SBarry Smith ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 2181aa26658SKarl Rupp } 2191aa26658SKarl Rupp if (fas->smoothd && fas->level != fas->levels - 1) { 2200298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 2211aa26658SKarl Rupp } 2221aa26658SKarl Rupp if (fas->smoothu && fas->level != fas->levels - 1) { 2230298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 2241aa26658SKarl Rupp } 22579d9a41aSPeter Brune } 22679d9a41aSPeter Brune 227534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 228534ebe21SPeter Brune if (fas->smoothd) { 229bc3f2f05SPeter Brune if (fas->level == 0 && fas->levels != 1) { 230365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); 231534ebe21SPeter Brune } else { 232365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 233534ebe21SPeter Brune } 2347fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr); 235534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); 2367601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2377601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr); 2386b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2396b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2406b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2416b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 242f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 2430dd27c6cSPeter Brune 2440dd27c6cSPeter Brune fas->smoothd->vec_sol = snes->vec_sol; 2450dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 2460dd27c6cSPeter Brune fas->smoothd->vec_sol_update = snes->vec_sol_update; 2470dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 2480dd27c6cSPeter Brune fas->smoothd->vec_func = snes->vec_func; 2490dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 2500dd27c6cSPeter Brune 2510dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2520dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr); 2530dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 254534ebe21SPeter Brune } 255534ebe21SPeter Brune 256534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 257534ebe21SPeter Brune if (fas->smoothu) { 258534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 259365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); 260534ebe21SPeter Brune } else { 261365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 262534ebe21SPeter Brune } 2637fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr); 264534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); 2657601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2667601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr); 2676b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2686b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2696b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2706b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 271f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 2720dd27c6cSPeter Brune 2730dd27c6cSPeter Brune fas->smoothu->vec_sol = snes->vec_sol; 2740dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 2750dd27c6cSPeter Brune fas->smoothu->vec_sol_update = snes->vec_sol_update; 2760dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 2770dd27c6cSPeter Brune fas->smoothu->vec_func = snes->vec_func; 2780dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 2790dd27c6cSPeter Brune 2800dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2810dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr); 2820dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2830dd27c6cSPeter Brune 284534ebe21SPeter Brune } 285d06165b7SPeter Brune 286ab8d36c9SPeter Brune if (next) { 28779d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 288ab8d36c9SPeter Brune if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 289ab8d36c9SPeter Brune if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 2907fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); 2917601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2927601faf0SJed Brown ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr); 2936b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2946b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2956b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2966b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 297f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 298ab8d36c9SPeter Brune ierr = SNESSetUp(next);CHKERRQ(ierr); 29979d9a41aSPeter Brune } 3006273346dSPeter Brune /* setup FAS work vectors */ 3016273346dSPeter Brune if (fas->galerkin) { 3026273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 3036273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 3046273346dSPeter Brune } 305421d9b32SPeter Brune PetscFunctionReturn(0); 306421d9b32SPeter Brune } 307421d9b32SPeter Brune 308421d9b32SPeter Brune #undef __FUNCT__ 309421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 310421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 311421d9b32SPeter Brune { 312ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 313ee78dd50SPeter Brune PetscInt levels = 1; 314*87f44e3fSPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE; 315421d9b32SPeter Brune PetscErrorCode ierr; 316ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 31707144faaSPeter Brune SNESFASType fastype; 318fde0ff24SPeter Brune const char *optionsprefix; 319f1c6b773SPeter Brune SNESLineSearch linesearch; 32066585501SPeter Brune PetscInt m, n_up, n_down; 321ab8d36c9SPeter Brune SNES next; 322ab8d36c9SPeter Brune PetscBool isFine; 323421d9b32SPeter Brune 324421d9b32SPeter Brune PetscFunctionBegin; 325ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 326c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 327ee78dd50SPeter Brune 328ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 329ab8d36c9SPeter Brune if (isFine) { 330ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 331c732cbdbSBarry Smith if (!flg && snes->dm) { 332c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 333c732cbdbSBarry Smith levels++; 334d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 335c732cbdbSBarry Smith } 3360298fd71SBarry Smith ierr = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr); 33707144faaSPeter Brune fastype = fas->fastype; 33807144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 33907144faaSPeter Brune if (flg) { 34007144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 34107144faaSPeter Brune } 342ee78dd50SPeter Brune 343fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 344ab8d36c9SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 345ab8d36c9SPeter Brune if (flg) { 346ab8d36c9SPeter Brune ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 347fde0ff24SPeter Brune } 348*87f44e3fSPeter Brune ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr); 349*87f44e3fSPeter Brune if (flg) { 350*87f44e3fSPeter Brune ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr); 351*87f44e3fSPeter Brune } 352fde0ff24SPeter Brune 353ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 354ab8d36c9SPeter Brune if (flg) { 355ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 356ab8d36c9SPeter Brune } 357ee78dd50SPeter Brune 358928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 359928e959bSPeter Brune ierr = PetscOptionsBool("-snes_fas_full_downsweep","Smooth on the initial upsweep for full FAS cycles","SNESFASFullSetDownSweep",fas->full_downsweep,&fas->full_downsweep,&flg);CHKERRQ(ierr); 360928e959bSPeter Brune if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);} 361928e959bSPeter Brune } 362928e959bSPeter Brune 36366585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 364162d76ddSPeter Brune 36566585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 366162d76ddSPeter Brune 367c8c899caSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 368c8c899caSPeter Brune if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 3690dd27c6cSPeter Brune 3700dd27c6cSPeter Brune flg = PETSC_FALSE; 3710dd27c6cSPeter Brune monflg = PETSC_TRUE; 3720dd27c6cSPeter Brune ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 3730dd27c6cSPeter Brune if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 374ab8d36c9SPeter Brune } 375ee78dd50SPeter Brune 376421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3778cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 378162d76ddSPeter Brune if (upflg) { 37966585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 380162d76ddSPeter Brune } 381162d76ddSPeter Brune if (downflg) { 38266585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 383162d76ddSPeter Brune } 384eff52c0eSPeter Brune 3859e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3869e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3879e764e56SPeter Brune if (!snes->linesearch) { 3887601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 3891a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 3909e764e56SPeter Brune } 3919e764e56SPeter Brune } 3929e764e56SPeter Brune 393ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 394ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 395ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 396421d9b32SPeter Brune PetscFunctionReturn(0); 397421d9b32SPeter Brune } 398421d9b32SPeter Brune 3999804daf3SBarry Smith #include <petscdraw.h> 400421d9b32SPeter Brune #undef __FUNCT__ 401421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 402421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 403421d9b32SPeter Brune { 404421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 405656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 406ab8d36c9SPeter Brune PetscInt i; 407421d9b32SPeter Brune PetscErrorCode ierr; 408ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 409421d9b32SPeter Brune 410421d9b32SPeter Brune PetscFunctionBegin; 411ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 412ab8d36c9SPeter Brune if (isFine) { 413251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 414656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 415421d9b32SPeter Brune if (iascii) { 416ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 417ab8d36c9SPeter Brune if (fas->galerkin) { 418ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 419421d9b32SPeter Brune } else { 420ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 421421d9b32SPeter Brune } 422ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 423ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 424ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 425ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 426ab8d36c9SPeter Brune if (!i) { 427ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 428421d9b32SPeter Brune } else { 429ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 430421d9b32SPeter Brune } 431ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 432166b3ea4SJed Brown if (smoothd) { 433ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 434166b3ea4SJed Brown } else { 435166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 436166b3ea4SJed Brown } 437ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 438ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 439ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 440ab8d36c9SPeter Brune } else if (i) { 441ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 442ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 443166b3ea4SJed Brown if (smoothu) { 444ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 445166b3ea4SJed Brown } else { 446166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 447166b3ea4SJed Brown } 448ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 449ab8d36c9SPeter Brune } 450ab8d36c9SPeter Brune } 451656ede7eSPeter Brune } else if (isdraw) { 452656ede7eSPeter Brune PetscDraw draw; 453b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 454656ede7eSPeter Brune SNES_FAS *curfas = fas; 455656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 456656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 457656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 458656ede7eSPeter Brune bottom = y - th; 459656ede7eSPeter Brune while (curfas) { 460b4375e8dSPeter Brune if (!curfas->smoothu) { 461656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 462656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 463656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 464b4375e8dSPeter Brune } else { 465b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 466b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 467b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 468b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 469b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 470b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 471b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 472b4375e8dSPeter Brune } 473656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 474656ede7eSPeter Brune bottom -= 5*th; 4751aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 4760298fd71SBarry Smith else curfas = NULL; 477656ede7eSPeter Brune } 478421d9b32SPeter Brune } 479ab8d36c9SPeter Brune } 480421d9b32SPeter Brune PetscFunctionReturn(0); 481421d9b32SPeter Brune } 482421d9b32SPeter Brune 483421d9b32SPeter Brune #undef __FUNCT__ 48491f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private" 48539bd7f45SPeter Brune /* 48639bd7f45SPeter Brune Defines the action of the downsmoother 48739bd7f45SPeter Brune */ 48891f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 489b9c2fdf1SPeter Brune { 49039bd7f45SPeter Brune PetscErrorCode ierr = 0; 491742fe5e2SPeter Brune SNESConvergedReason reason; 492ab8d36c9SPeter Brune Vec FPC; 493ab8d36c9SPeter Brune SNES smoothd; 4940dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 4956e111a19SKarl Rupp 496421d9b32SPeter Brune PetscFunctionBegin; 497ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 498e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 4990dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 500ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 5010dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 502742fe5e2SPeter Brune /* check convergence reason for the smoother */ 503ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 5041ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 505742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 506742fe5e2SPeter Brune PetscFunctionReturn(0); 507742fe5e2SPeter Brune } 5080298fd71SBarry Smith ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 5094b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 510b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr); 51139bd7f45SPeter Brune PetscFunctionReturn(0); 51239bd7f45SPeter Brune } 51339bd7f45SPeter Brune 51439bd7f45SPeter Brune 51539bd7f45SPeter Brune #undef __FUNCT__ 51691f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private" 51739bd7f45SPeter Brune /* 51807144faaSPeter Brune Defines the action of the upsmoother 51939bd7f45SPeter Brune */ 5200adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 5210adebc6cSBarry Smith { 52239bd7f45SPeter Brune PetscErrorCode ierr = 0; 52339bd7f45SPeter Brune SNESConvergedReason reason; 524ab8d36c9SPeter Brune Vec FPC; 525ab8d36c9SPeter Brune SNES smoothu; 5260dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 527ab8d36c9SPeter Brune 5286e111a19SKarl Rupp PetscFunctionBegin; 529ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 5300dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 531ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 5320dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 53339bd7f45SPeter Brune /* check convergence reason for the smoother */ 534ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 5351ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 53639bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 53739bd7f45SPeter Brune PetscFunctionReturn(0); 53839bd7f45SPeter Brune } 5390298fd71SBarry Smith ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 5404b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 541b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr); 54239bd7f45SPeter Brune PetscFunctionReturn(0); 54339bd7f45SPeter Brune } 54439bd7f45SPeter Brune 54539bd7f45SPeter Brune #undef __FUNCT__ 546938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 547938e4a01SJed Brown /*@ 548938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 549938e4a01SJed Brown 550938e4a01SJed Brown Collective 551938e4a01SJed Brown 552938e4a01SJed Brown Input Arguments: 553938e4a01SJed Brown . snes - SNESFAS 554938e4a01SJed Brown 555938e4a01SJed Brown Output Arguments: 556938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 557938e4a01SJed Brown 558938e4a01SJed Brown Level: developer 559938e4a01SJed Brown 560938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 561938e4a01SJed Brown @*/ 562938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 563938e4a01SJed Brown { 564938e4a01SJed Brown PetscErrorCode ierr; 565938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 566938e4a01SJed Brown 567938e4a01SJed Brown PetscFunctionBegin; 5681aa26658SKarl Rupp if (fas->rscale) { 5691aa26658SKarl Rupp ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 570f5af7f23SKarl Rupp } else if (fas->inject) { 5710298fd71SBarry Smith ierr = MatGetVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 572ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 573938e4a01SJed Brown PetscFunctionReturn(0); 574938e4a01SJed Brown } 575938e4a01SJed Brown 576e9923e8dSJed Brown #undef __FUNCT__ 577e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 578e9923e8dSJed Brown /*@ 579e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 580e9923e8dSJed Brown 581e9923e8dSJed Brown Collective 582e9923e8dSJed Brown 583e9923e8dSJed Brown Input Arguments: 584e9923e8dSJed Brown + fine - SNES from which to restrict 585e9923e8dSJed Brown - Xfine - vector to restrict 586e9923e8dSJed Brown 587e9923e8dSJed Brown Output Arguments: 588e9923e8dSJed Brown . Xcoarse - result of restriction 589e9923e8dSJed Brown 590e9923e8dSJed Brown Level: developer 591e9923e8dSJed Brown 592e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 593e9923e8dSJed Brown @*/ 594e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 595e9923e8dSJed Brown { 596e9923e8dSJed Brown PetscErrorCode ierr; 597e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 598e9923e8dSJed Brown 599e9923e8dSJed Brown PetscFunctionBegin; 600e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 601e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 602e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 603e9923e8dSJed Brown if (fas->inject) { 604e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 605e9923e8dSJed Brown } else { 606e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 607e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 608e9923e8dSJed Brown } 609e9923e8dSJed Brown PetscFunctionReturn(0); 610e9923e8dSJed Brown } 611e9923e8dSJed Brown 612e9923e8dSJed Brown #undef __FUNCT__ 6138c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection" 61439bd7f45SPeter Brune /* 61539bd7f45SPeter Brune 61639bd7f45SPeter Brune Performs the FAS coarse correction as: 61739bd7f45SPeter Brune 61839bd7f45SPeter Brune fine problem: F(x) = 0 61939bd7f45SPeter Brune coarse problem: F^c(x) = b^c 62039bd7f45SPeter Brune 62139bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 62239bd7f45SPeter Brune 62339bd7f45SPeter Brune */ 6240adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 6250adebc6cSBarry Smith { 62639bd7f45SPeter Brune PetscErrorCode ierr; 62739bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 62839bd7f45SPeter Brune SNESConvergedReason reason; 629ab8d36c9SPeter Brune SNES next; 630ab8d36c9SPeter Brune Mat restrct, interpolate; 6310dd27c6cSPeter Brune SNES_FAS *fasc; 6325fd66863SKarl Rupp 63339bd7f45SPeter Brune PetscFunctionBegin; 634ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 635ab8d36c9SPeter Brune if (next) { 6360dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 6370dd27c6cSPeter Brune 638ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 639ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 640ab8d36c9SPeter Brune 641ab8d36c9SPeter Brune X_c = next->vec_sol; 642ab8d36c9SPeter Brune Xo_c = next->work[0]; 643ab8d36c9SPeter Brune F_c = next->vec_func; 644ab8d36c9SPeter Brune B_c = next->vec_rhs; 645efe1f98aSPeter Brune 6460dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 647938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 648293a7e31SPeter Brune /* restrict the defect */ 649ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 6500dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 6510dd27c6cSPeter Brune 6520dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 653ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 6540dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6550dd27c6cSPeter Brune 6560dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 657e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 658b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 659e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 660ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 661ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 662c90fad12SPeter Brune 663c90fad12SPeter Brune /* recurse to the next level */ 664e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 665ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 666ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 667742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 668742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 669742fe5e2SPeter Brune PetscFunctionReturn(0); 670742fe5e2SPeter Brune } 671fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 672fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 6730dd27c6cSPeter Brune 6740dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 675ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 6760dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 677293a7e31SPeter Brune } 67839bd7f45SPeter Brune PetscFunctionReturn(0); 67939bd7f45SPeter Brune } 68039bd7f45SPeter Brune 68139bd7f45SPeter Brune #undef __FUNCT__ 6822cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive" 68339bd7f45SPeter Brune /* 68439bd7f45SPeter Brune 68539bd7f45SPeter Brune The additive cycle looks like: 68639bd7f45SPeter Brune 68707144faaSPeter Brune xhat = x 68807144faaSPeter Brune xhat = dS(x, b) 68907144faaSPeter Brune x = coarsecorrection(xhat, b_d) 69007144faaSPeter Brune x = x + nu*(xhat - x); 69139bd7f45SPeter Brune (optional) x = uS(x, b) 69239bd7f45SPeter Brune 69339bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 69439bd7f45SPeter Brune 69539bd7f45SPeter Brune */ 6960adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 6970adebc6cSBarry Smith { 69807144faaSPeter Brune Vec F, B, Xhat; 69922c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 70039bd7f45SPeter Brune PetscErrorCode ierr; 70107144faaSPeter Brune SNESConvergedReason reason; 70222c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 70322c1e704SPeter Brune PetscBool lssuccess; 704ab8d36c9SPeter Brune SNES next; 705ab8d36c9SPeter Brune Mat restrct, interpolate; 7060dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 7070dd27c6cSPeter Brune 70839bd7f45SPeter Brune PetscFunctionBegin; 709ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 71039bd7f45SPeter Brune F = snes->vec_func; 71139bd7f45SPeter Brune B = snes->vec_rhs; 712e7f468e7SPeter Brune Xhat = snes->work[1]; 71307144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 71407144faaSPeter Brune /* recurse first */ 715ab8d36c9SPeter Brune if (next) { 7160dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 717ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 718ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 7190dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 72007144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 7210dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 722c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 723ab8d36c9SPeter Brune X_c = next->vec_sol; 724ab8d36c9SPeter Brune Xo_c = next->work[0]; 725ab8d36c9SPeter Brune F_c = next->vec_func; 726ab8d36c9SPeter Brune B_c = next->vec_rhs; 72739bd7f45SPeter Brune 728938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 72907144faaSPeter Brune /* restrict the defect */ 730ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 73107144faaSPeter Brune 73207144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 7330dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 734ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 7350dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 736e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 737b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 738e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 73907144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 74007144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 74107144faaSPeter Brune 74207144faaSPeter Brune /* recurse */ 743e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 744ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 74507144faaSPeter Brune 74607144faaSPeter Brune /* smooth on this level */ 74791f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 74807144faaSPeter Brune 749ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 75007144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 75107144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 75207144faaSPeter Brune PetscFunctionReturn(0); 75307144faaSPeter Brune } 75407144faaSPeter Brune 75507144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 756c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 757ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 75807144faaSPeter Brune 759ddebd997SPeter Brune /* additive correction of the coarse direction*/ 760f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 761f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 7629e764e56SPeter Brune if (!lssuccess) { 7639e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 7649e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 7659e764e56SPeter Brune PetscFunctionReturn(0); 7669e764e56SPeter Brune } 7679e764e56SPeter Brune } 768b9c2fdf1SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 76907144faaSPeter Brune } else { 77091f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 77107144faaSPeter Brune } 77239bd7f45SPeter Brune PetscFunctionReturn(0); 77339bd7f45SPeter Brune } 77439bd7f45SPeter Brune 77539bd7f45SPeter Brune #undef __FUNCT__ 7762cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative" 77739bd7f45SPeter Brune /* 77839bd7f45SPeter Brune 77939bd7f45SPeter Brune Defines the FAS cycle as: 78039bd7f45SPeter Brune 78139bd7f45SPeter Brune fine problem: F(x) = 0 78239bd7f45SPeter Brune coarse problem: F^c(x) = b^c 78339bd7f45SPeter Brune 78439bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 78539bd7f45SPeter Brune 78639bd7f45SPeter Brune correction: 78739bd7f45SPeter Brune 78839bd7f45SPeter Brune x = x + I(x^c - Rx) 78939bd7f45SPeter Brune 79039bd7f45SPeter Brune */ 7910adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7920adebc6cSBarry Smith { 79339bd7f45SPeter Brune 79439bd7f45SPeter Brune PetscErrorCode ierr; 79539bd7f45SPeter Brune Vec F,B; 79634d65b3cSPeter Brune SNES next; 79739bd7f45SPeter Brune 79839bd7f45SPeter Brune PetscFunctionBegin; 79939bd7f45SPeter Brune F = snes->vec_func; 80039bd7f45SPeter Brune B = snes->vec_rhs; 80139bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 80234d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 80391f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 80434d65b3cSPeter Brune if (next) { 8058c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 80691f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 807fe6f9142SPeter Brune } 808fa9694d7SPeter Brune PetscFunctionReturn(0); 809421d9b32SPeter Brune } 810421d9b32SPeter Brune 811421d9b32SPeter Brune #undef __FUNCT__ 8128c94862eSPeter Brune #define __FUNCT__ "SNESFASCycleSetupPhase_Full" 8138c94862eSPeter Brune PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 8148c94862eSPeter Brune { 8158c94862eSPeter Brune SNES next; 8168c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8178c94862eSPeter Brune PetscBool isFine; 8188c94862eSPeter Brune PetscErrorCode ierr; 8198c94862eSPeter Brune 8208c94862eSPeter Brune PetscFunctionBegin; 8218c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 8228c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8238c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8248c94862eSPeter Brune fas->full_stage = 0; 8258c94862eSPeter Brune if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);} 8268c94862eSPeter Brune PetscFunctionReturn(0); 8278c94862eSPeter Brune } 8288c94862eSPeter Brune 8298c94862eSPeter Brune #undef __FUNCT__ 8308c94862eSPeter Brune #define __FUNCT__ "SNESFASCycle_Full" 8318c94862eSPeter Brune PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 8328c94862eSPeter Brune { 8338c94862eSPeter Brune PetscErrorCode ierr; 8348c94862eSPeter Brune Vec F,B; 8358c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8368c94862eSPeter Brune PetscBool isFine; 8378c94862eSPeter Brune SNES next; 8388c94862eSPeter Brune 8398c94862eSPeter Brune PetscFunctionBegin; 8408c94862eSPeter Brune F = snes->vec_func; 8418c94862eSPeter Brune B = snes->vec_rhs; 8428c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8438c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8448c94862eSPeter Brune 8458c94862eSPeter Brune if (isFine) { 8468c94862eSPeter Brune ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr); 8478c94862eSPeter Brune } 8488c94862eSPeter Brune 8498c94862eSPeter Brune if (fas->full_stage == 0) { 850928e959bSPeter Brune /* downsweep */ 8518c94862eSPeter Brune if (next) { 8528c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 853928e959bSPeter Brune if (fas->full_downsweep||isFine) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8548c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8558c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8568c94862eSPeter Brune if (fas->level != 1) next->max_its -= 1; 8578c94862eSPeter Brune } else { 8588c94862eSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8598c94862eSPeter Brune } 8608c94862eSPeter Brune fas->full_stage = 1; 8618c94862eSPeter Brune } else if (fas->full_stage == 1) { 8628c94862eSPeter Brune if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8638c94862eSPeter Brune if (next) { 8648c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8658c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8668c94862eSPeter Brune } 8678c94862eSPeter Brune } 8688c94862eSPeter Brune /* final v-cycle */ 8698c94862eSPeter Brune if (isFine) { 8708c94862eSPeter Brune if (next) { 8718c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8728c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8738c94862eSPeter Brune } 8748c94862eSPeter Brune } 8758c94862eSPeter Brune PetscFunctionReturn(0); 8768c94862eSPeter Brune } 8778c94862eSPeter Brune 8788c94862eSPeter Brune #undef __FUNCT__ 87934d65b3cSPeter Brune #define __FUNCT__ "SNESFASCycle_Kaskade" 88034d65b3cSPeter Brune PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 88134d65b3cSPeter Brune { 88234d65b3cSPeter Brune 88334d65b3cSPeter Brune PetscErrorCode ierr; 88434d65b3cSPeter Brune Vec F,B; 88534d65b3cSPeter Brune SNES next; 88634d65b3cSPeter Brune 88734d65b3cSPeter Brune PetscFunctionBegin; 88834d65b3cSPeter Brune F = snes->vec_func; 88934d65b3cSPeter Brune B = snes->vec_rhs; 89034d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 89134d65b3cSPeter Brune if (next) { 89234d65b3cSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 89334d65b3cSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 89434d65b3cSPeter Brune } else { 89534d65b3cSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 89634d65b3cSPeter Brune } 89734d65b3cSPeter Brune PetscFunctionReturn(0); 89834d65b3cSPeter Brune } 89934d65b3cSPeter Brune 90034d65b3cSPeter Brune #undef __FUNCT__ 901421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 902421d9b32SPeter Brune 903421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 904421d9b32SPeter Brune { 905fa9694d7SPeter Brune PetscErrorCode ierr; 906fe6f9142SPeter Brune PetscInt i, maxits; 907ddb5aff1SPeter Brune Vec X, F; 908fe6f9142SPeter Brune PetscReal fnorm; 909b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 910b17ce1afSJed Brown DM dm; 911e70c42e5SPeter Brune PetscBool isFine; 912b17ce1afSJed Brown 913421d9b32SPeter Brune PetscFunctionBegin; 914fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 915fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 916fa9694d7SPeter Brune X = snes->vec_sol; 917f5a6d4f9SBarry Smith F = snes->vec_func; 918293a7e31SPeter Brune 91918a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 920293a7e31SPeter Brune /*norm setup */ 921e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 922fe6f9142SPeter Brune snes->iter = 0; 923fe6f9142SPeter Brune snes->norm = 0.; 924e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 925e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 9260dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 927fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 9280dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 929fe6f9142SPeter Brune if (snes->domainerror) { 930fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 931fe6f9142SPeter Brune PetscFunctionReturn(0); 932fe6f9142SPeter Brune } 9331aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 934e4ed7901SPeter Brune 935fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 936189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 937189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 938189a9710SBarry Smith PetscFunctionReturn(0); 939189a9710SBarry Smith } 940e4ed7901SPeter Brune 941e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 942fe6f9142SPeter Brune snes->norm = fnorm; 943e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 944a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 945fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 946fe6f9142SPeter Brune 947fe6f9142SPeter Brune /* test convergence */ 948fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 949fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 950e4ed7901SPeter Brune 951b17ce1afSJed Brown 952b9c2fdf1SPeter Brune if (isFine) { 953b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 954b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 955b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 956b17ce1afSJed Brown DM dmcoarse; 957b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 958b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 959b17ce1afSJed Brown dm = dmcoarse; 960b17ce1afSJed Brown } 961b9c2fdf1SPeter Brune } 962b17ce1afSJed Brown 963fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 964fe6f9142SPeter Brune /* Call general purpose update function */ 965646217ecSPeter Brune 966fe6f9142SPeter Brune if (snes->ops->update) { 967fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 968fe6f9142SPeter Brune } 96907144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 97091f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 9718c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_ADDITIVE) { 97291f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 9738c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_FULL) { 9748c94862eSPeter Brune ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr); 97534d65b3cSPeter Brune } else if (fas->fastype ==SNES_FAS_KASKADE) { 97634d65b3cSPeter Brune ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr); 9778c94862eSPeter Brune } else { 9788c94862eSPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");CHKERRQ(ierr); 97907144faaSPeter Brune } 980742fe5e2SPeter Brune 981742fe5e2SPeter Brune /* check for FAS cycle divergence */ 9821aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 983b9c2fdf1SPeter Brune 984c90fad12SPeter Brune /* Monitor convergence */ 985e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 986c90fad12SPeter Brune snes->iter = i+1; 987e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 988a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 989c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 990c90fad12SPeter Brune /* Test for convergence */ 99166585501SPeter Brune if (isFine) { 992b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 993c90fad12SPeter Brune if (snes->reason) break; 994fe6f9142SPeter Brune } 99566585501SPeter Brune } 996fe6f9142SPeter Brune if (i == maxits) { 997fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 998fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 999fe6f9142SPeter Brune } 1000421d9b32SPeter Brune PetscFunctionReturn(0); 1001421d9b32SPeter Brune } 1002