1421d9b32SPeter Brune /* Defines the basic SNES object */ 2a055c8caSBarry Smith #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnes.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); 84416b707SBarry Smith extern PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,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 31d39b9438SMatthew G. Knepley . -snes_fas_full_downsweep<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 454f02bc6aSBarry Smith References: 4696a0c994SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 474f02bc6aSBarry Smith SIAM Review, 57(4), 2015 484f02bc6aSBarry Smith 49d3bc2379SPeter Brune .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 501fbfccc6SJed Brown M*/ 51421d9b32SPeter Brune 52421d9b32SPeter Brune #undef __FUNCT__ 53421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 548cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes) 55421d9b32SPeter Brune { 56421d9b32SPeter Brune SNES_FAS *fas; 57421d9b32SPeter Brune PetscErrorCode ierr; 58421d9b32SPeter Brune 59421d9b32SPeter Brune PetscFunctionBegin; 60421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 61421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 62421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 63421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 64421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 65421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 66421d9b32SPeter Brune 67ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 68ed020824SBarry Smith snes->usespc = PETSC_FALSE; 69ed020824SBarry Smith 7088976e71SPeter Brune if (!snes->tolerancesset) { 710e444f03SPeter Brune snes->max_funcs = 30000; 720e444f03SPeter Brune snes->max_its = 10000; 7388976e71SPeter Brune } 740e444f03SPeter Brune 75b00a9115SJed Brown ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr); 761aa26658SKarl Rupp 77421d9b32SPeter Brune snes->data = (void*) fas; 78421d9b32SPeter Brune fas->level = 0; 79293a7e31SPeter Brune fas->levels = 1; 80ee78dd50SPeter Brune fas->n_cycles = 1; 81ee78dd50SPeter Brune fas->max_up_it = 1; 82ee78dd50SPeter Brune fas->max_down_it = 1; 830298fd71SBarry Smith fas->smoothu = NULL; 840298fd71SBarry Smith fas->smoothd = NULL; 850298fd71SBarry Smith fas->next = NULL; 860298fd71SBarry Smith fas->previous = NULL; 87ab8d36c9SPeter Brune fas->fine = snes; 880298fd71SBarry Smith fas->interpolate = NULL; 890298fd71SBarry Smith fas->restrct = NULL; 900298fd71SBarry Smith fas->inject = NULL; 91cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 92ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 93928e959bSPeter Brune fas->full_downsweep = PETSC_FALSE; 940dd27c6cSPeter Brune 950298fd71SBarry Smith fas->eventsmoothsetup = 0; 960298fd71SBarry Smith fas->eventsmoothsolve = 0; 970298fd71SBarry Smith fas->eventresidual = 0; 980298fd71SBarry Smith fas->eventinterprestrict = 0; 99efe1f98aSPeter Brune PetscFunctionReturn(0); 100efe1f98aSPeter Brune } 101efe1f98aSPeter Brune 102421d9b32SPeter Brune #undef __FUNCT__ 103421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 104421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 105421d9b32SPeter Brune { 10677df8cc4SPeter Brune PetscErrorCode ierr = 0; 107421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 108421d9b32SPeter Brune 109421d9b32SPeter Brune PetscFunctionBegin; 110ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr); 111ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 1123dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 113bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 114bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 115bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 1167cd33a7bSPeter Brune if (fas->galerkin) { 1177cd33a7bSPeter Brune ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr); 1187cd33a7bSPeter Brune ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr); 1197cd33a7bSPeter Brune } 120d477d801SPeter Brune if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);} 121421d9b32SPeter Brune PetscFunctionReturn(0); 122421d9b32SPeter Brune } 123421d9b32SPeter Brune 124421d9b32SPeter Brune #undef __FUNCT__ 125421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 126421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 127421d9b32SPeter Brune { 128421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 129742fe5e2SPeter Brune PetscErrorCode ierr = 0; 130421d9b32SPeter Brune 131421d9b32SPeter Brune PetscFunctionBegin; 132421d9b32SPeter Brune /* recursively resets and then destroys */ 13379d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 1341aa26658SKarl Rupp if (fas->next) { 1351aa26658SKarl Rupp ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 1361aa26658SKarl Rupp } 137421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 138421d9b32SPeter Brune PetscFunctionReturn(0); 139421d9b32SPeter Brune } 140421d9b32SPeter Brune 141421d9b32SPeter Brune #undef __FUNCT__ 142421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 143421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 144421d9b32SPeter Brune { 14548bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 146421d9b32SPeter Brune PetscErrorCode ierr; 147d1adcc6fSPeter Brune PetscInt dm_levels; 1483dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 149ab8d36c9SPeter Brune SNES next; 150ab8d36c9SPeter Brune PetscBool isFine; 151f89ba88eSPeter Brune SNESLineSearch linesearch; 152f89ba88eSPeter Brune SNESLineSearch slinesearch; 153f89ba88eSPeter Brune void *lsprectx,*lspostctx; 1546b2b7091SBarry Smith PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 1556b2b7091SBarry Smith PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 156eff52c0eSPeter Brune 1576b2b7091SBarry Smith PetscFunctionBegin; 158ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 159ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 160d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 161d1adcc6fSPeter Brune dm_levels++; 162cc05f883SPeter Brune if (dm_levels > fas->levels) { 1632e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 1643dccd265SPeter Brune vec_sol = snes->vec_sol; 1653dccd265SPeter Brune vec_func = snes->vec_func; 1663dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 1673dccd265SPeter Brune vec_rhs = snes->vec_rhs; 1680298fd71SBarry Smith snes->vec_sol = NULL; 1690298fd71SBarry Smith snes->vec_func = NULL; 1700298fd71SBarry Smith snes->vec_sol_update = NULL; 1710298fd71SBarry Smith snes->vec_rhs = NULL; 1723dccd265SPeter Brune 1733dccd265SPeter Brune /* reset the number of levels */ 1740298fd71SBarry Smith ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr); 175cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 1763dccd265SPeter Brune 1773dccd265SPeter Brune snes->vec_sol = vec_sol; 1783dccd265SPeter Brune snes->vec_func = vec_func; 1793dccd265SPeter Brune snes->vec_rhs = vec_rhs; 1803dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 181d1adcc6fSPeter Brune } 182d1adcc6fSPeter Brune } 183ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 184ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 1853dccd265SPeter Brune 186fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 187cc05f883SPeter Brune 188ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 189ab8d36c9SPeter Brune if (!fas->smoothd) { 190ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 191ab8d36c9SPeter Brune } 192ab8d36c9SPeter Brune 19379d9a41aSPeter Brune if (snes->dm) { 194ab8d36c9SPeter Brune /* set the smoother DMs properly */ 195ab8d36c9SPeter Brune if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 196ab8d36c9SPeter Brune ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 19779d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 198ab8d36c9SPeter Brune if (next) { 19979d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 200ab8d36c9SPeter Brune if (!next->dm) { 201ce94432eSBarry Smith ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr); 202ab8d36c9SPeter Brune ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 20379d9a41aSPeter Brune } 20479d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 20579d9a41aSPeter Brune if (!fas->interpolate) { 206ab8d36c9SPeter Brune ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 207bccf9bb3SJed Brown if (!fas->restrct) { 208bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 20979d9a41aSPeter Brune fas->restrct = fas->interpolate; 21079d9a41aSPeter Brune } 211bccf9bb3SJed Brown } 21279d9a41aSPeter Brune /* set the injection from the DM */ 21379d9a41aSPeter Brune if (!fas->inject) { 2146dbf9973SLawrence Mitchell ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr); 21579d9a41aSPeter Brune } 21679d9a41aSPeter Brune } 21779d9a41aSPeter Brune } 21879d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 21979d9a41aSPeter Brune if (fas->galerkin) { 2201aa26658SKarl Rupp if (next) { 2210298fd71SBarry Smith ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 2221aa26658SKarl Rupp } 2231aa26658SKarl Rupp if (fas->smoothd && fas->level != fas->levels - 1) { 2240298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 2251aa26658SKarl Rupp } 2261aa26658SKarl Rupp if (fas->smoothu && fas->level != fas->levels - 1) { 2270298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 2281aa26658SKarl Rupp } 22979d9a41aSPeter Brune } 23079d9a41aSPeter Brune 231534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 232534ebe21SPeter Brune if (fas->smoothd) { 233bc3f2f05SPeter Brune if (fas->level == 0 && fas->levels != 1) { 234365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); 235534ebe21SPeter Brune } else { 236365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 237534ebe21SPeter Brune } 2387fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr); 239534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); 2407601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2417601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr); 2426b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2436b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2446b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2456b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 246f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 2470dd27c6cSPeter Brune 2480dd27c6cSPeter Brune fas->smoothd->vec_sol = snes->vec_sol; 2490dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 2500dd27c6cSPeter Brune fas->smoothd->vec_sol_update = snes->vec_sol_update; 2510dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 2520dd27c6cSPeter Brune fas->smoothd->vec_func = snes->vec_func; 2530dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 2540dd27c6cSPeter Brune 2550dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2560dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr); 2570dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 258534ebe21SPeter Brune } 259534ebe21SPeter Brune 260534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 261534ebe21SPeter Brune if (fas->smoothu) { 262534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 263365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); 264534ebe21SPeter Brune } else { 265365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 266534ebe21SPeter Brune } 2677fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr); 268534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); 2697601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2707601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr); 2716b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2726b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2736b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2746b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 275f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 2760dd27c6cSPeter Brune 2770dd27c6cSPeter Brune fas->smoothu->vec_sol = snes->vec_sol; 2780dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 2790dd27c6cSPeter Brune fas->smoothu->vec_sol_update = snes->vec_sol_update; 2800dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 2810dd27c6cSPeter Brune fas->smoothu->vec_func = snes->vec_func; 2820dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 2830dd27c6cSPeter Brune 2840dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2850dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr); 2860dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2870dd27c6cSPeter Brune 288534ebe21SPeter Brune } 289d06165b7SPeter Brune 290ab8d36c9SPeter Brune if (next) { 29179d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 292ab8d36c9SPeter Brune if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 293ab8d36c9SPeter Brune if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 2947fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); 2957601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2967601faf0SJed Brown ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr); 2976b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2986b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2996b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 3006b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 301f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 302ab8d36c9SPeter Brune ierr = SNESSetUp(next);CHKERRQ(ierr); 30379d9a41aSPeter Brune } 3046273346dSPeter Brune /* setup FAS work vectors */ 3056273346dSPeter Brune if (fas->galerkin) { 3066273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 3076273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 3086273346dSPeter Brune } 309421d9b32SPeter Brune PetscFunctionReturn(0); 310421d9b32SPeter Brune } 311421d9b32SPeter Brune 312421d9b32SPeter Brune #undef __FUNCT__ 313421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 3144416b707SBarry Smith PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes) 315421d9b32SPeter Brune { 316ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 317ee78dd50SPeter Brune PetscInt levels = 1; 31887f44e3fSPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE; 319421d9b32SPeter Brune PetscErrorCode ierr; 32007144faaSPeter Brune SNESFASType fastype; 321fde0ff24SPeter Brune const char *optionsprefix; 322f1c6b773SPeter Brune SNESLineSearch linesearch; 32366585501SPeter Brune PetscInt m, n_up, n_down; 324ab8d36c9SPeter Brune SNES next; 325ab8d36c9SPeter Brune PetscBool isFine; 326421d9b32SPeter Brune 327421d9b32SPeter Brune PetscFunctionBegin; 328ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 329e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNESFAS Options-----------------------------------");CHKERRQ(ierr); 330ee78dd50SPeter Brune 331ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 332ab8d36c9SPeter Brune if (isFine) { 333ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 334c732cbdbSBarry Smith if (!flg && snes->dm) { 335c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 336c732cbdbSBarry Smith levels++; 337d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 338c732cbdbSBarry Smith } 3390298fd71SBarry Smith ierr = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr); 34007144faaSPeter Brune fastype = fas->fastype; 34107144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 34207144faaSPeter Brune if (flg) { 34307144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 34407144faaSPeter Brune } 345ee78dd50SPeter Brune 346fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 347ab8d36c9SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 348ab8d36c9SPeter Brune if (flg) { 349ab8d36c9SPeter Brune ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 350fde0ff24SPeter Brune } 35187f44e3fSPeter Brune ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr); 35287f44e3fSPeter Brune if (flg) { 35387f44e3fSPeter Brune ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr); 35487f44e3fSPeter Brune } 355fde0ff24SPeter Brune 356ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 357ab8d36c9SPeter Brune if (flg) { 358ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 359ab8d36c9SPeter Brune } 360ee78dd50SPeter Brune 361928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 362928e959bSPeter 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); 363928e959bSPeter Brune if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);} 364928e959bSPeter Brune } 365928e959bSPeter Brune 36666585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 367162d76ddSPeter Brune 36866585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 369162d76ddSPeter Brune 370*d142ab34SLawrence Mitchell { 371*d142ab34SLawrence Mitchell PetscViewer viewer; 372*d142ab34SLawrence Mitchell PetscViewerFormat format; 373*d142ab34SLawrence Mitchell ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix, 374*d142ab34SLawrence Mitchell "-snes_fas_monitor",&viewer,&format,&monflg);CHKERRQ(ierr); 375*d142ab34SLawrence Mitchell if (monflg) { 376*d142ab34SLawrence Mitchell PetscViewerAndFormat *vf; 377*d142ab34SLawrence Mitchell ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 378*d142ab34SLawrence Mitchell ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 379*d142ab34SLawrence Mitchell ierr = SNESFASSetMonitor(snes,vf,PETSC_TRUE);CHKERRQ(ierr); 380*d142ab34SLawrence Mitchell } 381*d142ab34SLawrence Mitchell } 3820dd27c6cSPeter Brune flg = PETSC_FALSE; 3830dd27c6cSPeter Brune monflg = PETSC_TRUE; 3840dd27c6cSPeter Brune ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 3850dd27c6cSPeter Brune if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 386ab8d36c9SPeter Brune } 387ee78dd50SPeter Brune 388421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3898cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 390162d76ddSPeter Brune if (upflg) { 39166585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 392162d76ddSPeter Brune } 393162d76ddSPeter Brune if (downflg) { 39466585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 395162d76ddSPeter Brune } 396eff52c0eSPeter Brune 3979e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3989e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3999e764e56SPeter Brune if (!snes->linesearch) { 4007601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 4011a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 4029e764e56SPeter Brune } 4039e764e56SPeter Brune } 4049e764e56SPeter Brune 405ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 406ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 407ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 408421d9b32SPeter Brune PetscFunctionReturn(0); 409421d9b32SPeter Brune } 410421d9b32SPeter Brune 4119804daf3SBarry Smith #include <petscdraw.h> 412421d9b32SPeter Brune #undef __FUNCT__ 413421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 414421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 415421d9b32SPeter Brune { 416421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 417656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 418ab8d36c9SPeter Brune PetscInt i; 419421d9b32SPeter Brune PetscErrorCode ierr; 420ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 421421d9b32SPeter Brune 422421d9b32SPeter Brune PetscFunctionBegin; 423ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 424ab8d36c9SPeter Brune if (isFine) { 425251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 426656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 427421d9b32SPeter Brune if (iascii) { 428ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 429ab8d36c9SPeter Brune if (fas->galerkin) { 430ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 431421d9b32SPeter Brune } else { 432ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 433421d9b32SPeter Brune } 434ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 435ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 436ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 437ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 438ab8d36c9SPeter Brune if (!i) { 439ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 440421d9b32SPeter Brune } else { 441ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 442421d9b32SPeter Brune } 443ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 444166b3ea4SJed Brown if (smoothd) { 445ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 446166b3ea4SJed Brown } else { 447166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 448166b3ea4SJed Brown } 449ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 450ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 451ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 452ab8d36c9SPeter Brune } else if (i) { 453ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 454ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 455166b3ea4SJed Brown if (smoothu) { 456ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 457166b3ea4SJed Brown } else { 458166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 459166b3ea4SJed Brown } 460ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 461ab8d36c9SPeter Brune } 462ab8d36c9SPeter Brune } 463656ede7eSPeter Brune } else if (isdraw) { 464656ede7eSPeter Brune PetscDraw draw; 465b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 466656ede7eSPeter Brune SNES_FAS *curfas = fas; 467656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 468656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 469656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 470656ede7eSPeter Brune bottom = y - th; 471656ede7eSPeter Brune while (curfas) { 472b4375e8dSPeter Brune if (!curfas->smoothu) { 473656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 474656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 475656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 476b4375e8dSPeter Brune } else { 477b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 478b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 479b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 480b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 481b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 482b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 483b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 484b4375e8dSPeter Brune } 485656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 486656ede7eSPeter Brune bottom -= 5*th; 4871aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 4880298fd71SBarry Smith else curfas = NULL; 489656ede7eSPeter Brune } 490421d9b32SPeter Brune } 491ab8d36c9SPeter Brune } 492421d9b32SPeter Brune PetscFunctionReturn(0); 493421d9b32SPeter Brune } 494421d9b32SPeter Brune 495421d9b32SPeter Brune #undef __FUNCT__ 49691f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private" 49739bd7f45SPeter Brune /* 49839bd7f45SPeter Brune Defines the action of the downsmoother 49939bd7f45SPeter Brune */ 50091f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 501b9c2fdf1SPeter Brune { 50239bd7f45SPeter Brune PetscErrorCode ierr = 0; 503742fe5e2SPeter Brune SNESConvergedReason reason; 504ab8d36c9SPeter Brune Vec FPC; 505ab8d36c9SPeter Brune SNES smoothd; 5060dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 5076e111a19SKarl Rupp 508421d9b32SPeter Brune PetscFunctionBegin; 509ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 510e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 5110dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 512ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 5130dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 514742fe5e2SPeter Brune /* check convergence reason for the smoother */ 515ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 5161ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 517742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 518742fe5e2SPeter Brune PetscFunctionReturn(0); 519742fe5e2SPeter Brune } 5200298fd71SBarry Smith ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 5214b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 52271dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 52339bd7f45SPeter Brune PetscFunctionReturn(0); 52439bd7f45SPeter Brune } 52539bd7f45SPeter Brune 52639bd7f45SPeter Brune 52739bd7f45SPeter Brune #undef __FUNCT__ 52891f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private" 52939bd7f45SPeter Brune /* 53007144faaSPeter Brune Defines the action of the upsmoother 53139bd7f45SPeter Brune */ 5320adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 5330adebc6cSBarry Smith { 53439bd7f45SPeter Brune PetscErrorCode ierr = 0; 53539bd7f45SPeter Brune SNESConvergedReason reason; 536ab8d36c9SPeter Brune Vec FPC; 537ab8d36c9SPeter Brune SNES smoothu; 5380dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 539ab8d36c9SPeter Brune 5406e111a19SKarl Rupp PetscFunctionBegin; 541ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 5420dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 543ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 5440dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 54539bd7f45SPeter Brune /* check convergence reason for the smoother */ 546ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 5471ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 54839bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 54939bd7f45SPeter Brune PetscFunctionReturn(0); 55039bd7f45SPeter Brune } 5510298fd71SBarry Smith ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 5524b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 55371dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 55439bd7f45SPeter Brune PetscFunctionReturn(0); 55539bd7f45SPeter Brune } 55639bd7f45SPeter Brune 55739bd7f45SPeter Brune #undef __FUNCT__ 558938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 559938e4a01SJed Brown /*@ 560938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 561938e4a01SJed Brown 562938e4a01SJed Brown Collective 563938e4a01SJed Brown 564938e4a01SJed Brown Input Arguments: 565938e4a01SJed Brown . snes - SNESFAS 566938e4a01SJed Brown 567938e4a01SJed Brown Output Arguments: 568938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 569938e4a01SJed Brown 570938e4a01SJed Brown Level: developer 571938e4a01SJed Brown 572938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 573938e4a01SJed Brown @*/ 574938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 575938e4a01SJed Brown { 576938e4a01SJed Brown PetscErrorCode ierr; 577938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 578938e4a01SJed Brown 579938e4a01SJed Brown PetscFunctionBegin; 5801aa26658SKarl Rupp if (fas->rscale) { 5811aa26658SKarl Rupp ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 582f5af7f23SKarl Rupp } else if (fas->inject) { 5832a7a6963SBarry Smith ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 584ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 585938e4a01SJed Brown PetscFunctionReturn(0); 586938e4a01SJed Brown } 587938e4a01SJed Brown 588e9923e8dSJed Brown #undef __FUNCT__ 589e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 590e9923e8dSJed Brown /*@ 591e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 592e9923e8dSJed Brown 593e9923e8dSJed Brown Collective 594e9923e8dSJed Brown 595e9923e8dSJed Brown Input Arguments: 596e9923e8dSJed Brown + fine - SNES from which to restrict 597e9923e8dSJed Brown - Xfine - vector to restrict 598e9923e8dSJed Brown 599e9923e8dSJed Brown Output Arguments: 600e9923e8dSJed Brown . Xcoarse - result of restriction 601e9923e8dSJed Brown 602e9923e8dSJed Brown Level: developer 603e9923e8dSJed Brown 604e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 605e9923e8dSJed Brown @*/ 606e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 607e9923e8dSJed Brown { 608e9923e8dSJed Brown PetscErrorCode ierr; 609e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 610e9923e8dSJed Brown 611e9923e8dSJed Brown PetscFunctionBegin; 612e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 613e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 614e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 615e9923e8dSJed Brown if (fas->inject) { 616e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 617e9923e8dSJed Brown } else { 618e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 619e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 620e9923e8dSJed Brown } 621e9923e8dSJed Brown PetscFunctionReturn(0); 622e9923e8dSJed Brown } 623e9923e8dSJed Brown 624e9923e8dSJed Brown #undef __FUNCT__ 6258c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection" 62639bd7f45SPeter Brune /* 62739bd7f45SPeter Brune 62839bd7f45SPeter Brune Performs the FAS coarse correction as: 62939bd7f45SPeter Brune 630b5527d98SMatthew G. Knepley fine problem: F(x) = b 631b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c 63239bd7f45SPeter Brune 633b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 63439bd7f45SPeter Brune 63539bd7f45SPeter Brune */ 6360adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 6370adebc6cSBarry Smith { 63839bd7f45SPeter Brune PetscErrorCode ierr; 63939bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 64039bd7f45SPeter Brune SNESConvergedReason reason; 641ab8d36c9SPeter Brune SNES next; 642ab8d36c9SPeter Brune Mat restrct, interpolate; 6430dd27c6cSPeter Brune SNES_FAS *fasc; 6445fd66863SKarl Rupp 64539bd7f45SPeter Brune PetscFunctionBegin; 646ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 647ab8d36c9SPeter Brune if (next) { 6480dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 6490dd27c6cSPeter Brune 650ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 651ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 652ab8d36c9SPeter Brune 653ab8d36c9SPeter Brune X_c = next->vec_sol; 654ab8d36c9SPeter Brune Xo_c = next->work[0]; 655ab8d36c9SPeter Brune F_c = next->vec_func; 656ab8d36c9SPeter Brune B_c = next->vec_rhs; 657efe1f98aSPeter Brune 6580dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 659938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 6605609cbf2SMatthew G. Knepley /* restrict the defect: R(F(x) - b) */ 661ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 6620dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 6630dd27c6cSPeter Brune 6640dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6655609cbf2SMatthew G. Knepley /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */ 666ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 6670dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6680dd27c6cSPeter Brune 6690dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 670e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 671b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 672e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 673ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 674ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 675c90fad12SPeter Brune 676c90fad12SPeter Brune /* recurse to the next level */ 677e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 678ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 679ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 680742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 681742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 682742fe5e2SPeter Brune PetscFunctionReturn(0); 683742fe5e2SPeter Brune } 684fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 685fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 6860dd27c6cSPeter Brune 6870dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 688ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 6890dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 690293a7e31SPeter Brune } 69139bd7f45SPeter Brune PetscFunctionReturn(0); 69239bd7f45SPeter Brune } 69339bd7f45SPeter Brune 69439bd7f45SPeter Brune #undef __FUNCT__ 6952cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive" 69639bd7f45SPeter Brune /* 69739bd7f45SPeter Brune 69839bd7f45SPeter Brune The additive cycle looks like: 69939bd7f45SPeter Brune 70007144faaSPeter Brune xhat = x 70107144faaSPeter Brune xhat = dS(x, b) 70207144faaSPeter Brune x = coarsecorrection(xhat, b_d) 70307144faaSPeter Brune x = x + nu*(xhat - x); 70439bd7f45SPeter Brune (optional) x = uS(x, b) 70539bd7f45SPeter Brune 70639bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 70739bd7f45SPeter Brune 70839bd7f45SPeter Brune */ 7090adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 7100adebc6cSBarry Smith { 71107144faaSPeter Brune Vec F, B, Xhat; 71222c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 71339bd7f45SPeter Brune PetscErrorCode ierr; 71407144faaSPeter Brune SNESConvergedReason reason; 71522c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 716422a814eSBarry Smith SNESLineSearchReason lsresult; 717ab8d36c9SPeter Brune SNES next; 718ab8d36c9SPeter Brune Mat restrct, interpolate; 7190dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 7200dd27c6cSPeter Brune 72139bd7f45SPeter Brune PetscFunctionBegin; 722ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 72339bd7f45SPeter Brune F = snes->vec_func; 72439bd7f45SPeter Brune B = snes->vec_rhs; 725e7f468e7SPeter Brune Xhat = snes->work[1]; 72607144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 72707144faaSPeter Brune /* recurse first */ 728ab8d36c9SPeter Brune if (next) { 7290dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 730ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 731ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 7320dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 73307144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 7340dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 735c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 736ab8d36c9SPeter Brune X_c = next->vec_sol; 737ab8d36c9SPeter Brune Xo_c = next->work[0]; 738ab8d36c9SPeter Brune F_c = next->vec_func; 739ab8d36c9SPeter Brune B_c = next->vec_rhs; 74039bd7f45SPeter Brune 741938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 74207144faaSPeter Brune /* restrict the defect */ 743ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 74407144faaSPeter Brune 74507144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 7460dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 747ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 7480dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 749e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 750b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 751e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 75207144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 75307144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 75407144faaSPeter Brune 75507144faaSPeter Brune /* recurse */ 756e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 757ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 75807144faaSPeter Brune 75907144faaSPeter Brune /* smooth on this level */ 76091f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 76107144faaSPeter Brune 762ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 76307144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 76407144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 76507144faaSPeter Brune PetscFunctionReturn(0); 76607144faaSPeter Brune } 76707144faaSPeter Brune 76807144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 769c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 770ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 77107144faaSPeter Brune 772ddebd997SPeter Brune /* additive correction of the coarse direction*/ 773f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 774422a814eSBarry Smith ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr); 775422a814eSBarry Smith ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 776422a814eSBarry Smith if (lsresult) { 7779e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 7789e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 7799e764e56SPeter Brune PetscFunctionReturn(0); 7809e764e56SPeter Brune } 7819e764e56SPeter Brune } 78207144faaSPeter Brune } else { 78391f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 78407144faaSPeter Brune } 78539bd7f45SPeter Brune PetscFunctionReturn(0); 78639bd7f45SPeter Brune } 78739bd7f45SPeter Brune 78839bd7f45SPeter Brune #undef __FUNCT__ 7892cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative" 79039bd7f45SPeter Brune /* 79139bd7f45SPeter Brune 79239bd7f45SPeter Brune Defines the FAS cycle as: 79339bd7f45SPeter Brune 794b5527d98SMatthew G. Knepley fine problem: F(x) = b 79539bd7f45SPeter Brune coarse problem: F^c(x) = b^c 79639bd7f45SPeter Brune 797b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 79839bd7f45SPeter Brune 79939bd7f45SPeter Brune correction: 80039bd7f45SPeter Brune 80139bd7f45SPeter Brune x = x + I(x^c - Rx) 80239bd7f45SPeter Brune 80339bd7f45SPeter Brune */ 8040adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 8050adebc6cSBarry Smith { 80639bd7f45SPeter Brune 80739bd7f45SPeter Brune PetscErrorCode ierr; 80839bd7f45SPeter Brune Vec F,B; 80934d65b3cSPeter Brune SNES next; 81039bd7f45SPeter Brune 81139bd7f45SPeter Brune PetscFunctionBegin; 81239bd7f45SPeter Brune F = snes->vec_func; 81339bd7f45SPeter Brune B = snes->vec_rhs; 81439bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 81534d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 81691f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 81734d65b3cSPeter Brune if (next) { 8188c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 81991f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 820fe6f9142SPeter Brune } 821fa9694d7SPeter Brune PetscFunctionReturn(0); 822421d9b32SPeter Brune } 823421d9b32SPeter Brune 824421d9b32SPeter Brune #undef __FUNCT__ 8258c94862eSPeter Brune #define __FUNCT__ "SNESFASCycleSetupPhase_Full" 8268c94862eSPeter Brune PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 8278c94862eSPeter Brune { 8288c94862eSPeter Brune SNES next; 8298c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8308c94862eSPeter Brune PetscBool isFine; 8318c94862eSPeter Brune PetscErrorCode ierr; 8328c94862eSPeter Brune 8338c94862eSPeter Brune PetscFunctionBegin; 8348c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 8358c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8368c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8378c94862eSPeter Brune fas->full_stage = 0; 8388c94862eSPeter Brune if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);} 8398c94862eSPeter Brune PetscFunctionReturn(0); 8408c94862eSPeter Brune } 8418c94862eSPeter Brune 8428c94862eSPeter Brune #undef __FUNCT__ 8438c94862eSPeter Brune #define __FUNCT__ "SNESFASCycle_Full" 8448c94862eSPeter Brune PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 8458c94862eSPeter Brune { 8468c94862eSPeter Brune PetscErrorCode ierr; 8478c94862eSPeter Brune Vec F,B; 8488c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8498c94862eSPeter Brune PetscBool isFine; 8508c94862eSPeter Brune SNES next; 8518c94862eSPeter Brune 8528c94862eSPeter Brune PetscFunctionBegin; 8538c94862eSPeter Brune F = snes->vec_func; 8548c94862eSPeter Brune B = snes->vec_rhs; 8558c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8568c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8578c94862eSPeter Brune 8588c94862eSPeter Brune if (isFine) { 8598c94862eSPeter Brune ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr); 8608c94862eSPeter Brune } 8618c94862eSPeter Brune 8628c94862eSPeter Brune if (fas->full_stage == 0) { 863928e959bSPeter Brune /* downsweep */ 8648c94862eSPeter Brune if (next) { 8658c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 866928e959bSPeter Brune if (fas->full_downsweep||isFine) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 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 if (fas->level != 1) next->max_its -= 1; 8708c94862eSPeter Brune } else { 871a3a80b83SMatthew G. Knepley /* The smoother on the coarse level is the coarse solver */ 8728c94862eSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8738c94862eSPeter Brune } 8748c94862eSPeter Brune fas->full_stage = 1; 8758c94862eSPeter Brune } else if (fas->full_stage == 1) { 8768c94862eSPeter Brune if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8778c94862eSPeter Brune if (next) { 8788c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8798c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8808c94862eSPeter Brune } 8818c94862eSPeter Brune } 8828c94862eSPeter Brune /* final v-cycle */ 8838c94862eSPeter Brune if (isFine) { 8848c94862eSPeter Brune if (next) { 8858c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8868c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8878c94862eSPeter Brune } 8888c94862eSPeter Brune } 8898c94862eSPeter Brune PetscFunctionReturn(0); 8908c94862eSPeter Brune } 8918c94862eSPeter Brune 8928c94862eSPeter Brune #undef __FUNCT__ 89334d65b3cSPeter Brune #define __FUNCT__ "SNESFASCycle_Kaskade" 89434d65b3cSPeter Brune PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 89534d65b3cSPeter Brune { 89634d65b3cSPeter Brune PetscErrorCode ierr; 89734d65b3cSPeter Brune Vec F,B; 89834d65b3cSPeter Brune SNES next; 89934d65b3cSPeter Brune 90034d65b3cSPeter Brune PetscFunctionBegin; 90134d65b3cSPeter Brune F = snes->vec_func; 90234d65b3cSPeter Brune B = snes->vec_rhs; 90334d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 90434d65b3cSPeter Brune if (next) { 90534d65b3cSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 90634d65b3cSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 90734d65b3cSPeter Brune } else { 90834d65b3cSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 90934d65b3cSPeter Brune } 91034d65b3cSPeter Brune PetscFunctionReturn(0); 91134d65b3cSPeter Brune } 91234d65b3cSPeter Brune 913fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE; 914fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n" 915fffbeea8SBarry Smith " title = {Composing Scalable Nonlinear Algebraic Solvers},\n" 916fffbeea8SBarry Smith " author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n" 917fffbeea8SBarry Smith " year = 2013,\n" 918fffbeea8SBarry Smith " type = Preprint,\n" 919fffbeea8SBarry Smith " number = {ANL/MCS-P2010-0112},\n" 920fffbeea8SBarry Smith " institution = {Argonne National Laboratory}\n}\n"; 921fffbeea8SBarry Smith 92234d65b3cSPeter Brune #undef __FUNCT__ 923421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 924421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 925421d9b32SPeter Brune { 926fa9694d7SPeter Brune PetscErrorCode ierr; 927fe6f9142SPeter Brune PetscInt i, maxits; 928ddb5aff1SPeter Brune Vec X, F; 929fe6f9142SPeter Brune PetscReal fnorm; 930b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 931b17ce1afSJed Brown DM dm; 932e70c42e5SPeter Brune PetscBool isFine; 933b17ce1afSJed Brown 934421d9b32SPeter Brune PetscFunctionBegin; 935c579b300SPatrick Farrell 9366c4ed002SBarry Smith if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 937c579b300SPatrick Farrell 938fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 939fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 940fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 941fa9694d7SPeter Brune X = snes->vec_sol; 942f5a6d4f9SBarry Smith F = snes->vec_func; 943293a7e31SPeter Brune 94418a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 945293a7e31SPeter Brune /*norm setup */ 946e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 947fe6f9142SPeter Brune snes->iter = 0; 948fe6f9142SPeter Brune snes->norm = 0.; 949e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 950e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 9510dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 952fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 9530dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 9541aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 955e4ed7901SPeter Brune 956fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 957422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 958e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 959fe6f9142SPeter Brune snes->norm = fnorm; 960e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 961a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 962fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 963fe6f9142SPeter Brune 964fe6f9142SPeter Brune /* test convergence */ 965fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 966fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 967e4ed7901SPeter Brune 968b17ce1afSJed Brown 969b9c2fdf1SPeter Brune if (isFine) { 970b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 971b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 972b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 973b17ce1afSJed Brown DM dmcoarse; 974b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 975b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 976b17ce1afSJed Brown dm = dmcoarse; 977b17ce1afSJed Brown } 978b9c2fdf1SPeter Brune } 979b17ce1afSJed Brown 980fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 981fe6f9142SPeter Brune /* Call general purpose update function */ 982646217ecSPeter Brune 983fe6f9142SPeter Brune if (snes->ops->update) { 984fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 985fe6f9142SPeter Brune } 98607144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 98791f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 9888c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_ADDITIVE) { 98991f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 9908c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_FULL) { 9918c94862eSPeter Brune ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr); 99234d65b3cSPeter Brune } else if (fas->fastype ==SNES_FAS_KASKADE) { 99334d65b3cSPeter Brune ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr); 9946c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type"); 995742fe5e2SPeter Brune 996742fe5e2SPeter Brune /* check for FAS cycle divergence */ 9971aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 998b9c2fdf1SPeter Brune 999c90fad12SPeter Brune /* Monitor convergence */ 1000e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 1001c90fad12SPeter Brune snes->iter = i+1; 1002e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1003a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 1004c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1005c90fad12SPeter Brune /* Test for convergence */ 100666585501SPeter Brune if (isFine) { 1007b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1008c90fad12SPeter Brune if (snes->reason) break; 1009fe6f9142SPeter Brune } 101066585501SPeter Brune } 1011fe6f9142SPeter Brune if (i == maxits) { 1012fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1013fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1014fe6f9142SPeter Brune } 1015421d9b32SPeter Brune PetscFunctionReturn(0); 1016421d9b32SPeter Brune } 1017