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); 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 71b00a9115SJed 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); 1137cd33a7bSPeter Brune if (fas->galerkin) { 1147cd33a7bSPeter Brune ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr); 1157cd33a7bSPeter Brune ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr); 1167cd33a7bSPeter Brune } 117d477d801SPeter Brune if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);} 118421d9b32SPeter Brune PetscFunctionReturn(0); 119421d9b32SPeter Brune } 120421d9b32SPeter Brune 121421d9b32SPeter Brune #undef __FUNCT__ 122421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 123421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 124421d9b32SPeter Brune { 125421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 126742fe5e2SPeter Brune PetscErrorCode ierr = 0; 127421d9b32SPeter Brune 128421d9b32SPeter Brune PetscFunctionBegin; 129421d9b32SPeter Brune /* recursively resets and then destroys */ 13079d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 1311aa26658SKarl Rupp if (fas->next) { 1321aa26658SKarl Rupp ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 1331aa26658SKarl Rupp } 134421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 135421d9b32SPeter Brune PetscFunctionReturn(0); 136421d9b32SPeter Brune } 137421d9b32SPeter Brune 138421d9b32SPeter Brune #undef __FUNCT__ 139421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 140421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 141421d9b32SPeter Brune { 14248bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 143421d9b32SPeter Brune PetscErrorCode ierr; 144efe1f98aSPeter Brune VecScatter injscatter; 145d1adcc6fSPeter Brune PetscInt dm_levels; 1463dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 147ab8d36c9SPeter Brune SNES next; 148ab8d36c9SPeter Brune PetscBool isFine; 149f89ba88eSPeter Brune SNESLineSearch linesearch; 150f89ba88eSPeter Brune SNESLineSearch slinesearch; 151f89ba88eSPeter Brune void *lsprectx,*lspostctx; 1526b2b7091SBarry Smith PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 1536b2b7091SBarry Smith PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 154eff52c0eSPeter Brune 1556b2b7091SBarry Smith PetscFunctionBegin; 156ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 157ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 158d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 159d1adcc6fSPeter Brune dm_levels++; 160cc05f883SPeter Brune if (dm_levels > fas->levels) { 1612e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 1623dccd265SPeter Brune vec_sol = snes->vec_sol; 1633dccd265SPeter Brune vec_func = snes->vec_func; 1643dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 1653dccd265SPeter Brune vec_rhs = snes->vec_rhs; 1660298fd71SBarry Smith snes->vec_sol = NULL; 1670298fd71SBarry Smith snes->vec_func = NULL; 1680298fd71SBarry Smith snes->vec_sol_update = NULL; 1690298fd71SBarry Smith snes->vec_rhs = NULL; 1703dccd265SPeter Brune 1713dccd265SPeter Brune /* reset the number of levels */ 1720298fd71SBarry Smith ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr); 173cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 1743dccd265SPeter Brune 1753dccd265SPeter Brune snes->vec_sol = vec_sol; 1763dccd265SPeter Brune snes->vec_func = vec_func; 1773dccd265SPeter Brune snes->vec_rhs = vec_rhs; 1783dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 179d1adcc6fSPeter Brune } 180d1adcc6fSPeter Brune } 181ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 182ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 1833dccd265SPeter Brune 184fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 185cc05f883SPeter Brune 186ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 187ab8d36c9SPeter Brune if (!fas->smoothd) { 188ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 189ab8d36c9SPeter Brune } 190ab8d36c9SPeter Brune 19179d9a41aSPeter Brune if (snes->dm) { 192ab8d36c9SPeter Brune /* set the smoother DMs properly */ 193ab8d36c9SPeter Brune if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 194ab8d36c9SPeter Brune ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 19579d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 196ab8d36c9SPeter Brune if (next) { 19779d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 198ab8d36c9SPeter Brune if (!next->dm) { 199ce94432eSBarry Smith ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr); 200ab8d36c9SPeter Brune ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 20179d9a41aSPeter Brune } 20279d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 20379d9a41aSPeter Brune if (!fas->interpolate) { 204ab8d36c9SPeter Brune ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 205bccf9bb3SJed Brown if (!fas->restrct) { 206bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 20779d9a41aSPeter Brune fas->restrct = fas->interpolate; 20879d9a41aSPeter Brune } 209bccf9bb3SJed Brown } 21079d9a41aSPeter Brune /* set the injection from the DM */ 21179d9a41aSPeter Brune if (!fas->inject) { 212ab8d36c9SPeter Brune ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 213ce94432eSBarry Smith ierr = MatCreateScatter(PetscObjectComm((PetscObject)snes), injscatter, &fas->inject);CHKERRQ(ierr); 21479d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);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" 314421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(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; 320ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 32107144faaSPeter Brune SNESFASType fastype; 322fde0ff24SPeter Brune const char *optionsprefix; 323f1c6b773SPeter Brune SNESLineSearch linesearch; 32466585501SPeter Brune PetscInt m, n_up, n_down; 325ab8d36c9SPeter Brune SNES next; 326ab8d36c9SPeter Brune PetscBool isFine; 327421d9b32SPeter Brune 328421d9b32SPeter Brune PetscFunctionBegin; 329ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 330c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 331ee78dd50SPeter Brune 332ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 333ab8d36c9SPeter Brune if (isFine) { 334ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 335c732cbdbSBarry Smith if (!flg && snes->dm) { 336c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 337c732cbdbSBarry Smith levels++; 338d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 339c732cbdbSBarry Smith } 3400298fd71SBarry Smith ierr = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr); 34107144faaSPeter Brune fastype = fas->fastype; 34207144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 34307144faaSPeter Brune if (flg) { 34407144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 34507144faaSPeter Brune } 346ee78dd50SPeter Brune 347fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 348ab8d36c9SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 349ab8d36c9SPeter Brune if (flg) { 350ab8d36c9SPeter Brune ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 351fde0ff24SPeter Brune } 35287f44e3fSPeter Brune ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr); 35387f44e3fSPeter Brune if (flg) { 35487f44e3fSPeter Brune ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr); 35587f44e3fSPeter Brune } 356fde0ff24SPeter Brune 357ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 358ab8d36c9SPeter Brune if (flg) { 359ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 360ab8d36c9SPeter Brune } 361ee78dd50SPeter Brune 362928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 363928e959bSPeter 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); 364928e959bSPeter Brune if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);} 365928e959bSPeter Brune } 366928e959bSPeter Brune 36766585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 368162d76ddSPeter Brune 36966585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 370162d76ddSPeter Brune 371c8c899caSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 372c8c899caSPeter Brune if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 3730dd27c6cSPeter Brune 3740dd27c6cSPeter Brune flg = PETSC_FALSE; 3750dd27c6cSPeter Brune monflg = PETSC_TRUE; 3760dd27c6cSPeter Brune ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 3770dd27c6cSPeter Brune if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 378ab8d36c9SPeter Brune } 379ee78dd50SPeter Brune 380421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3818cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 382162d76ddSPeter Brune if (upflg) { 38366585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 384162d76ddSPeter Brune } 385162d76ddSPeter Brune if (downflg) { 38666585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 387162d76ddSPeter Brune } 388eff52c0eSPeter Brune 3899e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3909e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3919e764e56SPeter Brune if (!snes->linesearch) { 3927601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 3931a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 3949e764e56SPeter Brune } 3959e764e56SPeter Brune } 3969e764e56SPeter Brune 397ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 398ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 399ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 400421d9b32SPeter Brune PetscFunctionReturn(0); 401421d9b32SPeter Brune } 402421d9b32SPeter Brune 4039804daf3SBarry Smith #include <petscdraw.h> 404421d9b32SPeter Brune #undef __FUNCT__ 405421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 406421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 407421d9b32SPeter Brune { 408421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 409656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 410ab8d36c9SPeter Brune PetscInt i; 411421d9b32SPeter Brune PetscErrorCode ierr; 412ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 413421d9b32SPeter Brune 414421d9b32SPeter Brune PetscFunctionBegin; 415ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 416ab8d36c9SPeter Brune if (isFine) { 417251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 418656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 419421d9b32SPeter Brune if (iascii) { 420ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 421ab8d36c9SPeter Brune if (fas->galerkin) { 422ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 423421d9b32SPeter Brune } else { 424ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 425421d9b32SPeter Brune } 426ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 427ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 428ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 429ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 430ab8d36c9SPeter Brune if (!i) { 431ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 432421d9b32SPeter Brune } else { 433ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 434421d9b32SPeter Brune } 435ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 436166b3ea4SJed Brown if (smoothd) { 437ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 438166b3ea4SJed Brown } else { 439166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 440166b3ea4SJed Brown } 441ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 442ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 443ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 444ab8d36c9SPeter Brune } else if (i) { 445ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 446ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 447166b3ea4SJed Brown if (smoothu) { 448ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 449166b3ea4SJed Brown } else { 450166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 451166b3ea4SJed Brown } 452ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 453ab8d36c9SPeter Brune } 454ab8d36c9SPeter Brune } 455656ede7eSPeter Brune } else if (isdraw) { 456656ede7eSPeter Brune PetscDraw draw; 457b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 458656ede7eSPeter Brune SNES_FAS *curfas = fas; 459656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 460656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 461656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 462656ede7eSPeter Brune bottom = y - th; 463656ede7eSPeter Brune while (curfas) { 464b4375e8dSPeter Brune if (!curfas->smoothu) { 465656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 466656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 467656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 468b4375e8dSPeter Brune } else { 469b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 470b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 471b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 472b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 473b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 474b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 475b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 476b4375e8dSPeter Brune } 477656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 478656ede7eSPeter Brune bottom -= 5*th; 4791aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 4800298fd71SBarry Smith else curfas = NULL; 481656ede7eSPeter Brune } 482421d9b32SPeter Brune } 483ab8d36c9SPeter Brune } 484421d9b32SPeter Brune PetscFunctionReturn(0); 485421d9b32SPeter Brune } 486421d9b32SPeter Brune 487421d9b32SPeter Brune #undef __FUNCT__ 48891f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private" 48939bd7f45SPeter Brune /* 49039bd7f45SPeter Brune Defines the action of the downsmoother 49139bd7f45SPeter Brune */ 49291f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 493b9c2fdf1SPeter Brune { 49439bd7f45SPeter Brune PetscErrorCode ierr = 0; 495742fe5e2SPeter Brune SNESConvergedReason reason; 496ab8d36c9SPeter Brune Vec FPC; 497ab8d36c9SPeter Brune SNES smoothd; 4980dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 4996e111a19SKarl Rupp 500421d9b32SPeter Brune PetscFunctionBegin; 501ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 502e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 5030dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 504ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 5050dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 506742fe5e2SPeter Brune /* check convergence reason for the smoother */ 507ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 5081ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 509742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 510742fe5e2SPeter Brune PetscFunctionReturn(0); 511742fe5e2SPeter Brune } 5120298fd71SBarry Smith ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 5134b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 51471dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 51539bd7f45SPeter Brune PetscFunctionReturn(0); 51639bd7f45SPeter Brune } 51739bd7f45SPeter Brune 51839bd7f45SPeter Brune 51939bd7f45SPeter Brune #undef __FUNCT__ 52091f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private" 52139bd7f45SPeter Brune /* 52207144faaSPeter Brune Defines the action of the upsmoother 52339bd7f45SPeter Brune */ 5240adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 5250adebc6cSBarry Smith { 52639bd7f45SPeter Brune PetscErrorCode ierr = 0; 52739bd7f45SPeter Brune SNESConvergedReason reason; 528ab8d36c9SPeter Brune Vec FPC; 529ab8d36c9SPeter Brune SNES smoothu; 5300dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 531ab8d36c9SPeter Brune 5326e111a19SKarl Rupp PetscFunctionBegin; 533ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 5340dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 535ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 5360dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 53739bd7f45SPeter Brune /* check convergence reason for the smoother */ 538ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 5391ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 54039bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 54139bd7f45SPeter Brune PetscFunctionReturn(0); 54239bd7f45SPeter Brune } 5430298fd71SBarry Smith ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 5444b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 54571dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 54639bd7f45SPeter Brune PetscFunctionReturn(0); 54739bd7f45SPeter Brune } 54839bd7f45SPeter Brune 54939bd7f45SPeter Brune #undef __FUNCT__ 550938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 551938e4a01SJed Brown /*@ 552938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 553938e4a01SJed Brown 554938e4a01SJed Brown Collective 555938e4a01SJed Brown 556938e4a01SJed Brown Input Arguments: 557938e4a01SJed Brown . snes - SNESFAS 558938e4a01SJed Brown 559938e4a01SJed Brown Output Arguments: 560938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 561938e4a01SJed Brown 562938e4a01SJed Brown Level: developer 563938e4a01SJed Brown 564938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 565938e4a01SJed Brown @*/ 566938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 567938e4a01SJed Brown { 568938e4a01SJed Brown PetscErrorCode ierr; 569938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 570938e4a01SJed Brown 571938e4a01SJed Brown PetscFunctionBegin; 5721aa26658SKarl Rupp if (fas->rscale) { 5731aa26658SKarl Rupp ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 574f5af7f23SKarl Rupp } else if (fas->inject) { 5750298fd71SBarry Smith ierr = MatGetVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 576ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 577938e4a01SJed Brown PetscFunctionReturn(0); 578938e4a01SJed Brown } 579938e4a01SJed Brown 580e9923e8dSJed Brown #undef __FUNCT__ 581e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 582e9923e8dSJed Brown /*@ 583e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 584e9923e8dSJed Brown 585e9923e8dSJed Brown Collective 586e9923e8dSJed Brown 587e9923e8dSJed Brown Input Arguments: 588e9923e8dSJed Brown + fine - SNES from which to restrict 589e9923e8dSJed Brown - Xfine - vector to restrict 590e9923e8dSJed Brown 591e9923e8dSJed Brown Output Arguments: 592e9923e8dSJed Brown . Xcoarse - result of restriction 593e9923e8dSJed Brown 594e9923e8dSJed Brown Level: developer 595e9923e8dSJed Brown 596e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 597e9923e8dSJed Brown @*/ 598e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 599e9923e8dSJed Brown { 600e9923e8dSJed Brown PetscErrorCode ierr; 601e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 602e9923e8dSJed Brown 603e9923e8dSJed Brown PetscFunctionBegin; 604e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 605e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 606e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 607e9923e8dSJed Brown if (fas->inject) { 608e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 609e9923e8dSJed Brown } else { 610e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 611e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 612e9923e8dSJed Brown } 613e9923e8dSJed Brown PetscFunctionReturn(0); 614e9923e8dSJed Brown } 615e9923e8dSJed Brown 616e9923e8dSJed Brown #undef __FUNCT__ 6178c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection" 61839bd7f45SPeter Brune /* 61939bd7f45SPeter Brune 62039bd7f45SPeter Brune Performs the FAS coarse correction as: 62139bd7f45SPeter Brune 622b5527d98SMatthew G. Knepley fine problem: F(x) = b 623b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c 62439bd7f45SPeter Brune 625b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 62639bd7f45SPeter Brune 62739bd7f45SPeter Brune */ 6280adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 6290adebc6cSBarry Smith { 63039bd7f45SPeter Brune PetscErrorCode ierr; 63139bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 63239bd7f45SPeter Brune SNESConvergedReason reason; 633ab8d36c9SPeter Brune SNES next; 634ab8d36c9SPeter Brune Mat restrct, interpolate; 6350dd27c6cSPeter Brune SNES_FAS *fasc; 6365fd66863SKarl Rupp 63739bd7f45SPeter Brune PetscFunctionBegin; 638ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 639ab8d36c9SPeter Brune if (next) { 6400dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 6410dd27c6cSPeter Brune 642ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 643ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 644ab8d36c9SPeter Brune 645ab8d36c9SPeter Brune X_c = next->vec_sol; 646ab8d36c9SPeter Brune Xo_c = next->work[0]; 647ab8d36c9SPeter Brune F_c = next->vec_func; 648ab8d36c9SPeter Brune B_c = next->vec_rhs; 649efe1f98aSPeter Brune 6500dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 651938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 652293a7e31SPeter Brune /* restrict the defect */ 653ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 6540dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 6550dd27c6cSPeter Brune 6560dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 657ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 6580dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6590dd27c6cSPeter Brune 6600dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 661e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 662b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 663e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 664ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 665ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 666c90fad12SPeter Brune 667c90fad12SPeter Brune /* recurse to the next level */ 668e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 669ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 670ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 671742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 672742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 673742fe5e2SPeter Brune PetscFunctionReturn(0); 674742fe5e2SPeter Brune } 675fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 676fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 6770dd27c6cSPeter Brune 6780dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 679ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 6800dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 681293a7e31SPeter Brune } 68239bd7f45SPeter Brune PetscFunctionReturn(0); 68339bd7f45SPeter Brune } 68439bd7f45SPeter Brune 68539bd7f45SPeter Brune #undef __FUNCT__ 6862cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive" 68739bd7f45SPeter Brune /* 68839bd7f45SPeter Brune 68939bd7f45SPeter Brune The additive cycle looks like: 69039bd7f45SPeter Brune 69107144faaSPeter Brune xhat = x 69207144faaSPeter Brune xhat = dS(x, b) 69307144faaSPeter Brune x = coarsecorrection(xhat, b_d) 69407144faaSPeter Brune x = x + nu*(xhat - x); 69539bd7f45SPeter Brune (optional) x = uS(x, b) 69639bd7f45SPeter Brune 69739bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 69839bd7f45SPeter Brune 69939bd7f45SPeter Brune */ 7000adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 7010adebc6cSBarry Smith { 70207144faaSPeter Brune Vec F, B, Xhat; 70322c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 70439bd7f45SPeter Brune PetscErrorCode ierr; 70507144faaSPeter Brune SNESConvergedReason reason; 70622c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 70722c1e704SPeter Brune PetscBool lssuccess; 708ab8d36c9SPeter Brune SNES next; 709ab8d36c9SPeter Brune Mat restrct, interpolate; 7100dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 7110dd27c6cSPeter Brune 71239bd7f45SPeter Brune PetscFunctionBegin; 713ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 71439bd7f45SPeter Brune F = snes->vec_func; 71539bd7f45SPeter Brune B = snes->vec_rhs; 716e7f468e7SPeter Brune Xhat = snes->work[1]; 71707144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 71807144faaSPeter Brune /* recurse first */ 719ab8d36c9SPeter Brune if (next) { 7200dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 721ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 722ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 7230dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 72407144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 7250dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 726c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 727ab8d36c9SPeter Brune X_c = next->vec_sol; 728ab8d36c9SPeter Brune Xo_c = next->work[0]; 729ab8d36c9SPeter Brune F_c = next->vec_func; 730ab8d36c9SPeter Brune B_c = next->vec_rhs; 73139bd7f45SPeter Brune 732938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 73307144faaSPeter Brune /* restrict the defect */ 734ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 73507144faaSPeter Brune 73607144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 7370dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 738ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 7390dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 740e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 741b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 742e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 74307144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 74407144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 74507144faaSPeter Brune 74607144faaSPeter Brune /* recurse */ 747e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 748ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 74907144faaSPeter Brune 75007144faaSPeter Brune /* smooth on this level */ 75191f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 75207144faaSPeter Brune 753ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 75407144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 75507144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 75607144faaSPeter Brune PetscFunctionReturn(0); 75707144faaSPeter Brune } 75807144faaSPeter Brune 75907144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 760c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 761ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 76207144faaSPeter Brune 763ddebd997SPeter Brune /* additive correction of the coarse direction*/ 764f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 765f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 7669e764e56SPeter Brune if (!lssuccess) { 7679e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 7689e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 7699e764e56SPeter Brune PetscFunctionReturn(0); 7709e764e56SPeter Brune } 7719e764e56SPeter Brune } 772b9c2fdf1SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 77307144faaSPeter Brune } else { 77491f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 77507144faaSPeter Brune } 77639bd7f45SPeter Brune PetscFunctionReturn(0); 77739bd7f45SPeter Brune } 77839bd7f45SPeter Brune 77939bd7f45SPeter Brune #undef __FUNCT__ 7802cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative" 78139bd7f45SPeter Brune /* 78239bd7f45SPeter Brune 78339bd7f45SPeter Brune Defines the FAS cycle as: 78439bd7f45SPeter Brune 785b5527d98SMatthew G. Knepley fine problem: F(x) = b 78639bd7f45SPeter Brune coarse problem: F^c(x) = b^c 78739bd7f45SPeter Brune 788b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 78939bd7f45SPeter Brune 79039bd7f45SPeter Brune correction: 79139bd7f45SPeter Brune 79239bd7f45SPeter Brune x = x + I(x^c - Rx) 79339bd7f45SPeter Brune 79439bd7f45SPeter Brune */ 7950adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7960adebc6cSBarry Smith { 79739bd7f45SPeter Brune 79839bd7f45SPeter Brune PetscErrorCode ierr; 79939bd7f45SPeter Brune Vec F,B; 80034d65b3cSPeter Brune SNES next; 80139bd7f45SPeter Brune 80239bd7f45SPeter Brune PetscFunctionBegin; 80339bd7f45SPeter Brune F = snes->vec_func; 80439bd7f45SPeter Brune B = snes->vec_rhs; 80539bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 80634d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 80791f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 80834d65b3cSPeter Brune if (next) { 8098c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 81091f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 811fe6f9142SPeter Brune } 812fa9694d7SPeter Brune PetscFunctionReturn(0); 813421d9b32SPeter Brune } 814421d9b32SPeter Brune 815421d9b32SPeter Brune #undef __FUNCT__ 8168c94862eSPeter Brune #define __FUNCT__ "SNESFASCycleSetupPhase_Full" 8178c94862eSPeter Brune PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 8188c94862eSPeter Brune { 8198c94862eSPeter Brune SNES next; 8208c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8218c94862eSPeter Brune PetscBool isFine; 8228c94862eSPeter Brune PetscErrorCode ierr; 8238c94862eSPeter Brune 8248c94862eSPeter Brune PetscFunctionBegin; 8258c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 8268c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8278c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8288c94862eSPeter Brune fas->full_stage = 0; 8298c94862eSPeter Brune if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);} 8308c94862eSPeter Brune PetscFunctionReturn(0); 8318c94862eSPeter Brune } 8328c94862eSPeter Brune 8338c94862eSPeter Brune #undef __FUNCT__ 8348c94862eSPeter Brune #define __FUNCT__ "SNESFASCycle_Full" 8358c94862eSPeter Brune PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 8368c94862eSPeter Brune { 8378c94862eSPeter Brune PetscErrorCode ierr; 8388c94862eSPeter Brune Vec F,B; 8398c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8408c94862eSPeter Brune PetscBool isFine; 8418c94862eSPeter Brune SNES next; 8428c94862eSPeter Brune 8438c94862eSPeter Brune PetscFunctionBegin; 8448c94862eSPeter Brune F = snes->vec_func; 8458c94862eSPeter Brune B = snes->vec_rhs; 8468c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8478c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8488c94862eSPeter Brune 8498c94862eSPeter Brune if (isFine) { 8508c94862eSPeter Brune ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr); 8518c94862eSPeter Brune } 8528c94862eSPeter Brune 8538c94862eSPeter Brune if (fas->full_stage == 0) { 854928e959bSPeter Brune /* downsweep */ 8558c94862eSPeter Brune if (next) { 8568c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 857928e959bSPeter Brune if (fas->full_downsweep||isFine) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8588c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8598c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8608c94862eSPeter Brune if (fas->level != 1) next->max_its -= 1; 8618c94862eSPeter Brune } else { 8628c94862eSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8638c94862eSPeter Brune } 8648c94862eSPeter Brune fas->full_stage = 1; 8658c94862eSPeter Brune } else if (fas->full_stage == 1) { 8668c94862eSPeter Brune if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8678c94862eSPeter Brune if (next) { 8688c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8698c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8708c94862eSPeter Brune } 8718c94862eSPeter Brune } 8728c94862eSPeter Brune /* final v-cycle */ 8738c94862eSPeter Brune if (isFine) { 8748c94862eSPeter Brune if (next) { 8758c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8768c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8778c94862eSPeter Brune } 8788c94862eSPeter Brune } 8798c94862eSPeter Brune PetscFunctionReturn(0); 8808c94862eSPeter Brune } 8818c94862eSPeter Brune 8828c94862eSPeter Brune #undef __FUNCT__ 88334d65b3cSPeter Brune #define __FUNCT__ "SNESFASCycle_Kaskade" 88434d65b3cSPeter Brune PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 88534d65b3cSPeter Brune { 88634d65b3cSPeter Brune PetscErrorCode ierr; 88734d65b3cSPeter Brune Vec F,B; 88834d65b3cSPeter Brune SNES next; 88934d65b3cSPeter Brune 89034d65b3cSPeter Brune PetscFunctionBegin; 89134d65b3cSPeter Brune F = snes->vec_func; 89234d65b3cSPeter Brune B = snes->vec_rhs; 89334d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 89434d65b3cSPeter Brune if (next) { 89534d65b3cSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 89634d65b3cSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 89734d65b3cSPeter Brune } else { 89834d65b3cSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 89934d65b3cSPeter Brune } 90034d65b3cSPeter Brune PetscFunctionReturn(0); 90134d65b3cSPeter Brune } 90234d65b3cSPeter Brune 903fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE; 904fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n" 905fffbeea8SBarry Smith " title = {Composing Scalable Nonlinear Algebraic Solvers},\n" 906fffbeea8SBarry Smith " author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n" 907fffbeea8SBarry Smith " year = 2013,\n" 908fffbeea8SBarry Smith " type = Preprint,\n" 909fffbeea8SBarry Smith " number = {ANL/MCS-P2010-0112},\n" 910fffbeea8SBarry Smith " institution = {Argonne National Laboratory}\n}\n"; 911fffbeea8SBarry Smith 91234d65b3cSPeter Brune #undef __FUNCT__ 913421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 914421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 915421d9b32SPeter Brune { 916fa9694d7SPeter Brune PetscErrorCode ierr; 917fe6f9142SPeter Brune PetscInt i, maxits; 918ddb5aff1SPeter Brune Vec X, F; 919fe6f9142SPeter Brune PetscReal fnorm; 920b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 921b17ce1afSJed Brown DM dm; 922e70c42e5SPeter Brune PetscBool isFine; 923b17ce1afSJed Brown 924421d9b32SPeter Brune PetscFunctionBegin; 925fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 926fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 927fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 928fa9694d7SPeter Brune X = snes->vec_sol; 929f5a6d4f9SBarry Smith F = snes->vec_func; 930293a7e31SPeter Brune 93118a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 932293a7e31SPeter Brune /*norm setup */ 933e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 934fe6f9142SPeter Brune snes->iter = 0; 935fe6f9142SPeter Brune snes->norm = 0.; 936e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 937e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 9380dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 939fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 9400dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 941fe6f9142SPeter Brune if (snes->domainerror) { 942fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 943fe6f9142SPeter Brune PetscFunctionReturn(0); 944fe6f9142SPeter Brune } 9451aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 946e4ed7901SPeter Brune 947fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 948189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 949189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 950189a9710SBarry Smith PetscFunctionReturn(0); 951189a9710SBarry Smith } 952e4ed7901SPeter Brune 953e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 954fe6f9142SPeter Brune snes->norm = fnorm; 955e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 956a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 957fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 958fe6f9142SPeter Brune 959fe6f9142SPeter Brune /* test convergence */ 960fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 961fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 962e4ed7901SPeter Brune 963b17ce1afSJed Brown 964b9c2fdf1SPeter Brune if (isFine) { 965b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 966b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 967b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 968b17ce1afSJed Brown DM dmcoarse; 969b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 970b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 971b17ce1afSJed Brown dm = dmcoarse; 972b17ce1afSJed Brown } 973b9c2fdf1SPeter Brune } 974b17ce1afSJed Brown 975fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 976fe6f9142SPeter Brune /* Call general purpose update function */ 977646217ecSPeter Brune 978fe6f9142SPeter Brune if (snes->ops->update) { 979fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 980fe6f9142SPeter Brune } 98107144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 98291f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 9838c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_ADDITIVE) { 98491f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 9858c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_FULL) { 9868c94862eSPeter Brune ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr); 98734d65b3cSPeter Brune } else if (fas->fastype ==SNES_FAS_KASKADE) { 98834d65b3cSPeter Brune ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr); 9898c94862eSPeter Brune } else { 990*4e619c20SJed Brown SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type"); 99107144faaSPeter Brune } 992742fe5e2SPeter Brune 993742fe5e2SPeter Brune /* check for FAS cycle divergence */ 9941aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 995b9c2fdf1SPeter Brune 996c90fad12SPeter Brune /* Monitor convergence */ 997e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 998c90fad12SPeter Brune snes->iter = i+1; 999e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 1000a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 1001c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1002c90fad12SPeter Brune /* Test for convergence */ 100366585501SPeter Brune if (isFine) { 1004b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1005c90fad12SPeter Brune if (snes->reason) break; 1006fe6f9142SPeter Brune } 100766585501SPeter Brune } 1008fe6f9142SPeter Brune if (i == maxits) { 1009fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1010fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1011fe6f9142SPeter Brune } 1012421d9b32SPeter Brune PetscFunctionReturn(0); 1013421d9b32SPeter Brune } 1014