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 71*b00a9115SJed Brown ierr = PetscNewLog(snes,&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; 3144d26bfa5SPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = 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 } 348fde0ff24SPeter Brune 349ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 350ab8d36c9SPeter Brune if (flg) { 351ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 352ab8d36c9SPeter Brune } 353ee78dd50SPeter Brune 354928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 355928e959bSPeter 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); 356928e959bSPeter Brune if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);} 357928e959bSPeter Brune } 358928e959bSPeter Brune 35966585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 360162d76ddSPeter Brune 36166585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 362162d76ddSPeter Brune 363c8c899caSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 364c8c899caSPeter Brune if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 3650dd27c6cSPeter Brune 3660dd27c6cSPeter Brune flg = PETSC_FALSE; 3670dd27c6cSPeter Brune monflg = PETSC_TRUE; 3680dd27c6cSPeter Brune ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 3690dd27c6cSPeter Brune if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 370ab8d36c9SPeter Brune } 371ee78dd50SPeter Brune 372421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3738cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 374162d76ddSPeter Brune if (upflg) { 37566585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 376162d76ddSPeter Brune } 377162d76ddSPeter Brune if (downflg) { 37866585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 379162d76ddSPeter Brune } 380eff52c0eSPeter Brune 3819e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3829e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3839e764e56SPeter Brune if (!snes->linesearch) { 3847601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 3851a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 3869e764e56SPeter Brune } 3879e764e56SPeter Brune } 3889e764e56SPeter Brune 389ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 390ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 391ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 392421d9b32SPeter Brune PetscFunctionReturn(0); 393421d9b32SPeter Brune } 394421d9b32SPeter Brune 3959804daf3SBarry Smith #include <petscdraw.h> 396421d9b32SPeter Brune #undef __FUNCT__ 397421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 398421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 399421d9b32SPeter Brune { 400421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 401656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 402ab8d36c9SPeter Brune PetscInt i; 403421d9b32SPeter Brune PetscErrorCode ierr; 404ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 405421d9b32SPeter Brune 406421d9b32SPeter Brune PetscFunctionBegin; 407ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 408ab8d36c9SPeter Brune if (isFine) { 409251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 410656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 411421d9b32SPeter Brune if (iascii) { 412ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 413ab8d36c9SPeter Brune if (fas->galerkin) { 414ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 415421d9b32SPeter Brune } else { 416ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 417421d9b32SPeter Brune } 418ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 419ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 420ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 421ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 422ab8d36c9SPeter Brune if (!i) { 423ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 424421d9b32SPeter Brune } else { 425ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 426421d9b32SPeter Brune } 427ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 428166b3ea4SJed Brown if (smoothd) { 429ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 430166b3ea4SJed Brown } else { 431166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 432166b3ea4SJed Brown } 433ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 434ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 435ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 436ab8d36c9SPeter Brune } else if (i) { 437ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 438ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 439166b3ea4SJed Brown if (smoothu) { 440ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 441166b3ea4SJed Brown } else { 442166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 443166b3ea4SJed Brown } 444ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 445ab8d36c9SPeter Brune } 446ab8d36c9SPeter Brune } 447656ede7eSPeter Brune } else if (isdraw) { 448656ede7eSPeter Brune PetscDraw draw; 449b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 450656ede7eSPeter Brune SNES_FAS *curfas = fas; 451656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 452656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 453656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 454656ede7eSPeter Brune bottom = y - th; 455656ede7eSPeter Brune while (curfas) { 456b4375e8dSPeter Brune if (!curfas->smoothu) { 457656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 458656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 459656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 460b4375e8dSPeter Brune } else { 461b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 462b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 463b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 464b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 465b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 466b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 467b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 468b4375e8dSPeter Brune } 469656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 470656ede7eSPeter Brune bottom -= 5*th; 4711aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 4720298fd71SBarry Smith else curfas = NULL; 473656ede7eSPeter Brune } 474421d9b32SPeter Brune } 475ab8d36c9SPeter Brune } 476421d9b32SPeter Brune PetscFunctionReturn(0); 477421d9b32SPeter Brune } 478421d9b32SPeter Brune 479421d9b32SPeter Brune #undef __FUNCT__ 48091f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private" 48139bd7f45SPeter Brune /* 48239bd7f45SPeter Brune Defines the action of the downsmoother 48339bd7f45SPeter Brune */ 48491f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 485b9c2fdf1SPeter Brune { 48639bd7f45SPeter Brune PetscErrorCode ierr = 0; 487742fe5e2SPeter Brune SNESConvergedReason reason; 488ab8d36c9SPeter Brune Vec FPC; 489ab8d36c9SPeter Brune SNES smoothd; 4900dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 4916e111a19SKarl Rupp 492421d9b32SPeter Brune PetscFunctionBegin; 493ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 494e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 4950dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 496ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 4970dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 498742fe5e2SPeter Brune /* check convergence reason for the smoother */ 499ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 5001ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 501742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 502742fe5e2SPeter Brune PetscFunctionReturn(0); 503742fe5e2SPeter Brune } 5040298fd71SBarry Smith ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 5054b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 506b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr); 50739bd7f45SPeter Brune PetscFunctionReturn(0); 50839bd7f45SPeter Brune } 50939bd7f45SPeter Brune 51039bd7f45SPeter Brune 51139bd7f45SPeter Brune #undef __FUNCT__ 51291f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private" 51339bd7f45SPeter Brune /* 51407144faaSPeter Brune Defines the action of the upsmoother 51539bd7f45SPeter Brune */ 5160adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 5170adebc6cSBarry Smith { 51839bd7f45SPeter Brune PetscErrorCode ierr = 0; 51939bd7f45SPeter Brune SNESConvergedReason reason; 520ab8d36c9SPeter Brune Vec FPC; 521ab8d36c9SPeter Brune SNES smoothu; 5220dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 523ab8d36c9SPeter Brune 5246e111a19SKarl Rupp PetscFunctionBegin; 525ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 5260dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 527ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 5280dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 52939bd7f45SPeter Brune /* check convergence reason for the smoother */ 530ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 5311ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 53239bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 53339bd7f45SPeter Brune PetscFunctionReturn(0); 53439bd7f45SPeter Brune } 5350298fd71SBarry Smith ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 5364b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 537b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr); 53839bd7f45SPeter Brune PetscFunctionReturn(0); 53939bd7f45SPeter Brune } 54039bd7f45SPeter Brune 54139bd7f45SPeter Brune #undef __FUNCT__ 542938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 543938e4a01SJed Brown /*@ 544938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 545938e4a01SJed Brown 546938e4a01SJed Brown Collective 547938e4a01SJed Brown 548938e4a01SJed Brown Input Arguments: 549938e4a01SJed Brown . snes - SNESFAS 550938e4a01SJed Brown 551938e4a01SJed Brown Output Arguments: 552938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 553938e4a01SJed Brown 554938e4a01SJed Brown Level: developer 555938e4a01SJed Brown 556938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 557938e4a01SJed Brown @*/ 558938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 559938e4a01SJed Brown { 560938e4a01SJed Brown PetscErrorCode ierr; 561938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 562938e4a01SJed Brown 563938e4a01SJed Brown PetscFunctionBegin; 5641aa26658SKarl Rupp if (fas->rscale) { 5651aa26658SKarl Rupp ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 566f5af7f23SKarl Rupp } else if (fas->inject) { 5670298fd71SBarry Smith ierr = MatGetVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 568ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 569938e4a01SJed Brown PetscFunctionReturn(0); 570938e4a01SJed Brown } 571938e4a01SJed Brown 572e9923e8dSJed Brown #undef __FUNCT__ 573e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 574e9923e8dSJed Brown /*@ 575e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 576e9923e8dSJed Brown 577e9923e8dSJed Brown Collective 578e9923e8dSJed Brown 579e9923e8dSJed Brown Input Arguments: 580e9923e8dSJed Brown + fine - SNES from which to restrict 581e9923e8dSJed Brown - Xfine - vector to restrict 582e9923e8dSJed Brown 583e9923e8dSJed Brown Output Arguments: 584e9923e8dSJed Brown . Xcoarse - result of restriction 585e9923e8dSJed Brown 586e9923e8dSJed Brown Level: developer 587e9923e8dSJed Brown 588e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 589e9923e8dSJed Brown @*/ 590e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 591e9923e8dSJed Brown { 592e9923e8dSJed Brown PetscErrorCode ierr; 593e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 594e9923e8dSJed Brown 595e9923e8dSJed Brown PetscFunctionBegin; 596e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 597e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 598e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 599e9923e8dSJed Brown if (fas->inject) { 600e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 601e9923e8dSJed Brown } else { 602e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 603e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 604e9923e8dSJed Brown } 605e9923e8dSJed Brown PetscFunctionReturn(0); 606e9923e8dSJed Brown } 607e9923e8dSJed Brown 608e9923e8dSJed Brown #undef __FUNCT__ 6098c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection" 61039bd7f45SPeter Brune /* 61139bd7f45SPeter Brune 61239bd7f45SPeter Brune Performs the FAS coarse correction as: 61339bd7f45SPeter Brune 61439bd7f45SPeter Brune fine problem: F(x) = 0 61539bd7f45SPeter Brune coarse problem: F^c(x) = b^c 61639bd7f45SPeter Brune 61739bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 61839bd7f45SPeter Brune 61939bd7f45SPeter Brune */ 6200adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 6210adebc6cSBarry Smith { 62239bd7f45SPeter Brune PetscErrorCode ierr; 62339bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 62439bd7f45SPeter Brune SNESConvergedReason reason; 625ab8d36c9SPeter Brune SNES next; 626ab8d36c9SPeter Brune Mat restrct, interpolate; 6270dd27c6cSPeter Brune SNES_FAS *fasc; 6285fd66863SKarl Rupp 62939bd7f45SPeter Brune PetscFunctionBegin; 630ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 631ab8d36c9SPeter Brune if (next) { 6320dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 6330dd27c6cSPeter Brune 634ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 635ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 636ab8d36c9SPeter Brune 637ab8d36c9SPeter Brune X_c = next->vec_sol; 638ab8d36c9SPeter Brune Xo_c = next->work[0]; 639ab8d36c9SPeter Brune F_c = next->vec_func; 640ab8d36c9SPeter Brune B_c = next->vec_rhs; 641efe1f98aSPeter Brune 6420dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 643938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 644293a7e31SPeter Brune /* restrict the defect */ 645ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 6460dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 6470dd27c6cSPeter Brune 6480dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 649ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 6500dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6510dd27c6cSPeter Brune 6520dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 653e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 654b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 655e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 656ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 657ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 658c90fad12SPeter Brune 659c90fad12SPeter Brune /* recurse to the next level */ 660e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 661ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 662ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 663742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 664742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 665742fe5e2SPeter Brune PetscFunctionReturn(0); 666742fe5e2SPeter Brune } 667fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 668fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 6690dd27c6cSPeter Brune 6700dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 671ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 6720dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 673293a7e31SPeter Brune } 67439bd7f45SPeter Brune PetscFunctionReturn(0); 67539bd7f45SPeter Brune } 67639bd7f45SPeter Brune 67739bd7f45SPeter Brune #undef __FUNCT__ 6782cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive" 67939bd7f45SPeter Brune /* 68039bd7f45SPeter Brune 68139bd7f45SPeter Brune The additive cycle looks like: 68239bd7f45SPeter Brune 68307144faaSPeter Brune xhat = x 68407144faaSPeter Brune xhat = dS(x, b) 68507144faaSPeter Brune x = coarsecorrection(xhat, b_d) 68607144faaSPeter Brune x = x + nu*(xhat - x); 68739bd7f45SPeter Brune (optional) x = uS(x, b) 68839bd7f45SPeter Brune 68939bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 69039bd7f45SPeter Brune 69139bd7f45SPeter Brune */ 6920adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 6930adebc6cSBarry Smith { 69407144faaSPeter Brune Vec F, B, Xhat; 69522c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 69639bd7f45SPeter Brune PetscErrorCode ierr; 69707144faaSPeter Brune SNESConvergedReason reason; 69822c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 69922c1e704SPeter Brune PetscBool lssuccess; 700ab8d36c9SPeter Brune SNES next; 701ab8d36c9SPeter Brune Mat restrct, interpolate; 7020dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 7030dd27c6cSPeter Brune 70439bd7f45SPeter Brune PetscFunctionBegin; 705ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 70639bd7f45SPeter Brune F = snes->vec_func; 70739bd7f45SPeter Brune B = snes->vec_rhs; 708e7f468e7SPeter Brune Xhat = snes->work[1]; 70907144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 71007144faaSPeter Brune /* recurse first */ 711ab8d36c9SPeter Brune if (next) { 7120dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 713ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 714ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 7150dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 71607144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 7170dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 718c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 719ab8d36c9SPeter Brune X_c = next->vec_sol; 720ab8d36c9SPeter Brune Xo_c = next->work[0]; 721ab8d36c9SPeter Brune F_c = next->vec_func; 722ab8d36c9SPeter Brune B_c = next->vec_rhs; 72339bd7f45SPeter Brune 724938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 72507144faaSPeter Brune /* restrict the defect */ 726ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 72707144faaSPeter Brune 72807144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 7290dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 730ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 7310dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 732e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 733b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 734e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 73507144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 73607144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 73707144faaSPeter Brune 73807144faaSPeter Brune /* recurse */ 739e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 740ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 74107144faaSPeter Brune 74207144faaSPeter Brune /* smooth on this level */ 74391f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 74407144faaSPeter Brune 745ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 74607144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 74707144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 74807144faaSPeter Brune PetscFunctionReturn(0); 74907144faaSPeter Brune } 75007144faaSPeter Brune 75107144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 752c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 753ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 75407144faaSPeter Brune 755ddebd997SPeter Brune /* additive correction of the coarse direction*/ 756f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 757f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 7589e764e56SPeter Brune if (!lssuccess) { 7599e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 7609e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 7619e764e56SPeter Brune PetscFunctionReturn(0); 7629e764e56SPeter Brune } 7639e764e56SPeter Brune } 764b9c2fdf1SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 76507144faaSPeter Brune } else { 76691f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 76707144faaSPeter Brune } 76839bd7f45SPeter Brune PetscFunctionReturn(0); 76939bd7f45SPeter Brune } 77039bd7f45SPeter Brune 77139bd7f45SPeter Brune #undef __FUNCT__ 7722cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative" 77339bd7f45SPeter Brune /* 77439bd7f45SPeter Brune 77539bd7f45SPeter Brune Defines the FAS cycle as: 77639bd7f45SPeter Brune 77739bd7f45SPeter Brune fine problem: F(x) = 0 77839bd7f45SPeter Brune coarse problem: F^c(x) = b^c 77939bd7f45SPeter Brune 78039bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 78139bd7f45SPeter Brune 78239bd7f45SPeter Brune correction: 78339bd7f45SPeter Brune 78439bd7f45SPeter Brune x = x + I(x^c - Rx) 78539bd7f45SPeter Brune 78639bd7f45SPeter Brune */ 7870adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7880adebc6cSBarry Smith { 78939bd7f45SPeter Brune 79039bd7f45SPeter Brune PetscErrorCode ierr; 79139bd7f45SPeter Brune Vec F,B; 79234d65b3cSPeter Brune SNES next; 79339bd7f45SPeter Brune 79439bd7f45SPeter Brune PetscFunctionBegin; 79539bd7f45SPeter Brune F = snes->vec_func; 79639bd7f45SPeter Brune B = snes->vec_rhs; 79739bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 79834d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 79991f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 80034d65b3cSPeter Brune if (next) { 8018c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 80291f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 803fe6f9142SPeter Brune } 804fa9694d7SPeter Brune PetscFunctionReturn(0); 805421d9b32SPeter Brune } 806421d9b32SPeter Brune 807421d9b32SPeter Brune #undef __FUNCT__ 8088c94862eSPeter Brune #define __FUNCT__ "SNESFASCycleSetupPhase_Full" 8098c94862eSPeter Brune PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 8108c94862eSPeter Brune { 8118c94862eSPeter Brune SNES next; 8128c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8138c94862eSPeter Brune PetscBool isFine; 8148c94862eSPeter Brune PetscErrorCode ierr; 8158c94862eSPeter Brune 8168c94862eSPeter Brune PetscFunctionBegin; 8178c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 8188c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8198c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8208c94862eSPeter Brune fas->full_stage = 0; 8218c94862eSPeter Brune if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);} 8228c94862eSPeter Brune PetscFunctionReturn(0); 8238c94862eSPeter Brune } 8248c94862eSPeter Brune 8258c94862eSPeter Brune #undef __FUNCT__ 8268c94862eSPeter Brune #define __FUNCT__ "SNESFASCycle_Full" 8278c94862eSPeter Brune PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 8288c94862eSPeter Brune { 8298c94862eSPeter Brune PetscErrorCode ierr; 8308c94862eSPeter Brune Vec F,B; 8318c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8328c94862eSPeter Brune PetscBool isFine; 8338c94862eSPeter Brune SNES next; 8348c94862eSPeter Brune 8358c94862eSPeter Brune PetscFunctionBegin; 8368c94862eSPeter Brune F = snes->vec_func; 8378c94862eSPeter Brune B = snes->vec_rhs; 8388c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8398c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8408c94862eSPeter Brune 8418c94862eSPeter Brune if (isFine) { 8428c94862eSPeter Brune ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr); 8438c94862eSPeter Brune } 8448c94862eSPeter Brune 8458c94862eSPeter Brune if (fas->full_stage == 0) { 846928e959bSPeter Brune /* downsweep */ 8478c94862eSPeter Brune if (next) { 8488c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 849928e959bSPeter Brune if (fas->full_downsweep||isFine) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8508c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8518c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8528c94862eSPeter Brune if (fas->level != 1) next->max_its -= 1; 8538c94862eSPeter Brune } else { 8548c94862eSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8558c94862eSPeter Brune } 8568c94862eSPeter Brune fas->full_stage = 1; 8578c94862eSPeter Brune } else if (fas->full_stage == 1) { 8588c94862eSPeter Brune if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8598c94862eSPeter Brune if (next) { 8608c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8618c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8628c94862eSPeter Brune } 8638c94862eSPeter Brune } 8648c94862eSPeter Brune /* final v-cycle */ 8658c94862eSPeter Brune if (isFine) { 8668c94862eSPeter Brune if (next) { 8678c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8688c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8698c94862eSPeter Brune } 8708c94862eSPeter Brune } 8718c94862eSPeter Brune PetscFunctionReturn(0); 8728c94862eSPeter Brune } 8738c94862eSPeter Brune 8748c94862eSPeter Brune #undef __FUNCT__ 87534d65b3cSPeter Brune #define __FUNCT__ "SNESFASCycle_Kaskade" 87634d65b3cSPeter Brune PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 87734d65b3cSPeter Brune { 87834d65b3cSPeter Brune 87934d65b3cSPeter Brune PetscErrorCode ierr; 88034d65b3cSPeter Brune Vec F,B; 88134d65b3cSPeter Brune SNES next; 88234d65b3cSPeter Brune 88334d65b3cSPeter Brune PetscFunctionBegin; 88434d65b3cSPeter Brune F = snes->vec_func; 88534d65b3cSPeter Brune B = snes->vec_rhs; 88634d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 88734d65b3cSPeter Brune if (next) { 88834d65b3cSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 88934d65b3cSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 89034d65b3cSPeter Brune } else { 89134d65b3cSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 89234d65b3cSPeter Brune } 89334d65b3cSPeter Brune PetscFunctionReturn(0); 89434d65b3cSPeter Brune } 89534d65b3cSPeter Brune 89634d65b3cSPeter Brune #undef __FUNCT__ 897421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 898421d9b32SPeter Brune 899421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 900421d9b32SPeter Brune { 901fa9694d7SPeter Brune PetscErrorCode ierr; 902fe6f9142SPeter Brune PetscInt i, maxits; 903ddb5aff1SPeter Brune Vec X, F; 904fe6f9142SPeter Brune PetscReal fnorm; 905b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 906b17ce1afSJed Brown DM dm; 907e70c42e5SPeter Brune PetscBool isFine; 908b17ce1afSJed Brown 909421d9b32SPeter Brune PetscFunctionBegin; 910fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 911fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 912fa9694d7SPeter Brune X = snes->vec_sol; 913f5a6d4f9SBarry Smith F = snes->vec_func; 914293a7e31SPeter Brune 91518a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 916293a7e31SPeter Brune /*norm setup */ 917e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 918fe6f9142SPeter Brune snes->iter = 0; 919fe6f9142SPeter Brune snes->norm = 0.; 920e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 921e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 9220dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 923fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 9240dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 925fe6f9142SPeter Brune if (snes->domainerror) { 926fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 927fe6f9142SPeter Brune PetscFunctionReturn(0); 928fe6f9142SPeter Brune } 9291aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 930e4ed7901SPeter Brune 931fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 932189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 933189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 934189a9710SBarry Smith PetscFunctionReturn(0); 935189a9710SBarry Smith } 936e4ed7901SPeter Brune 937e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 938fe6f9142SPeter Brune snes->norm = fnorm; 939e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 940a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 941fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 942fe6f9142SPeter Brune 943fe6f9142SPeter Brune /* test convergence */ 944fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 945fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 946e4ed7901SPeter Brune 947b17ce1afSJed Brown 948b9c2fdf1SPeter Brune if (isFine) { 949b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 950b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 951b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 952b17ce1afSJed Brown DM dmcoarse; 953b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 954b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 955b17ce1afSJed Brown dm = dmcoarse; 956b17ce1afSJed Brown } 957b9c2fdf1SPeter Brune } 958b17ce1afSJed Brown 959fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 960fe6f9142SPeter Brune /* Call general purpose update function */ 961646217ecSPeter Brune 962fe6f9142SPeter Brune if (snes->ops->update) { 963fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 964fe6f9142SPeter Brune } 96507144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 96691f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 9678c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_ADDITIVE) { 96891f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 9698c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_FULL) { 9708c94862eSPeter Brune ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr); 97134d65b3cSPeter Brune } else if (fas->fastype ==SNES_FAS_KASKADE) { 97234d65b3cSPeter Brune ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr); 9738c94862eSPeter Brune } else { 9748c94862eSPeter Brune SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");CHKERRQ(ierr); 97507144faaSPeter Brune } 976742fe5e2SPeter Brune 977742fe5e2SPeter Brune /* check for FAS cycle divergence */ 9781aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 979b9c2fdf1SPeter Brune 980c90fad12SPeter Brune /* Monitor convergence */ 981e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 982c90fad12SPeter Brune snes->iter = i+1; 983e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 984a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 985c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 986c90fad12SPeter Brune /* Test for convergence */ 98766585501SPeter Brune if (isFine) { 988b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 989c90fad12SPeter Brune if (snes->reason) break; 990fe6f9142SPeter Brune } 99166585501SPeter Brune } 992fe6f9142SPeter Brune if (i == maxits) { 993fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 994fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 995fe6f9142SPeter Brune } 996421d9b32SPeter Brune PetscFunctionReturn(0); 997421d9b32SPeter Brune } 998