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) { 212*6dbf9973SLawrence Mitchell ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr); 21379d9a41aSPeter Brune } 21479d9a41aSPeter Brune } 21579d9a41aSPeter Brune } 21679d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 21779d9a41aSPeter Brune if (fas->galerkin) { 2181aa26658SKarl Rupp if (next) { 2190298fd71SBarry Smith ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 2201aa26658SKarl Rupp } 2211aa26658SKarl Rupp if (fas->smoothd && fas->level != fas->levels - 1) { 2220298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 2231aa26658SKarl Rupp } 2241aa26658SKarl Rupp if (fas->smoothu && fas->level != fas->levels - 1) { 2250298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 2261aa26658SKarl Rupp } 22779d9a41aSPeter Brune } 22879d9a41aSPeter Brune 229534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 230534ebe21SPeter Brune if (fas->smoothd) { 231bc3f2f05SPeter Brune if (fas->level == 0 && fas->levels != 1) { 232365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); 233534ebe21SPeter Brune } else { 234365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 235534ebe21SPeter Brune } 2367fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr); 237534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); 2387601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2397601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr); 2406b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2416b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2426b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2436b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 244f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 2450dd27c6cSPeter Brune 2460dd27c6cSPeter Brune fas->smoothd->vec_sol = snes->vec_sol; 2470dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 2480dd27c6cSPeter Brune fas->smoothd->vec_sol_update = snes->vec_sol_update; 2490dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 2500dd27c6cSPeter Brune fas->smoothd->vec_func = snes->vec_func; 2510dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 2520dd27c6cSPeter Brune 2530dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2540dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr); 2550dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 256534ebe21SPeter Brune } 257534ebe21SPeter Brune 258534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 259534ebe21SPeter Brune if (fas->smoothu) { 260534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 261365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); 262534ebe21SPeter Brune } else { 263365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 264534ebe21SPeter Brune } 2657fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr); 266534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); 2677601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2687601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr); 2696b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2706b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2716b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2726b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 273f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 2740dd27c6cSPeter Brune 2750dd27c6cSPeter Brune fas->smoothu->vec_sol = snes->vec_sol; 2760dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 2770dd27c6cSPeter Brune fas->smoothu->vec_sol_update = snes->vec_sol_update; 2780dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 2790dd27c6cSPeter Brune fas->smoothu->vec_func = snes->vec_func; 2800dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 2810dd27c6cSPeter Brune 2820dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2830dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr); 2840dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2850dd27c6cSPeter Brune 286534ebe21SPeter Brune } 287d06165b7SPeter Brune 288ab8d36c9SPeter Brune if (next) { 28979d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 290ab8d36c9SPeter Brune if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 291ab8d36c9SPeter Brune if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 2927fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); 2937601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2947601faf0SJed Brown ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr); 2956b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2966b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2976b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2986b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 299f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 300ab8d36c9SPeter Brune ierr = SNESSetUp(next);CHKERRQ(ierr); 30179d9a41aSPeter Brune } 3026273346dSPeter Brune /* setup FAS work vectors */ 3036273346dSPeter Brune if (fas->galerkin) { 3046273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 3056273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 3066273346dSPeter Brune } 307421d9b32SPeter Brune PetscFunctionReturn(0); 308421d9b32SPeter Brune } 309421d9b32SPeter Brune 310421d9b32SPeter Brune #undef __FUNCT__ 311421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 312421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes) 313421d9b32SPeter Brune { 314ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 315ee78dd50SPeter Brune PetscInt levels = 1; 31687f44e3fSPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE; 317421d9b32SPeter Brune PetscErrorCode ierr; 318ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 31907144faaSPeter Brune SNESFASType fastype; 320fde0ff24SPeter Brune const char *optionsprefix; 321f1c6b773SPeter Brune SNESLineSearch linesearch; 32266585501SPeter Brune PetscInt m, n_up, n_down; 323ab8d36c9SPeter Brune SNES next; 324ab8d36c9SPeter Brune PetscBool isFine; 325421d9b32SPeter Brune 326421d9b32SPeter Brune PetscFunctionBegin; 327ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 328c90fad12SPeter Brune ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr); 329ee78dd50SPeter Brune 330ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 331ab8d36c9SPeter Brune if (isFine) { 332ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 333c732cbdbSBarry Smith if (!flg && snes->dm) { 334c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 335c732cbdbSBarry Smith levels++; 336d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 337c732cbdbSBarry Smith } 3380298fd71SBarry Smith ierr = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr); 33907144faaSPeter Brune fastype = fas->fastype; 34007144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 34107144faaSPeter Brune if (flg) { 34207144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 34307144faaSPeter Brune } 344ee78dd50SPeter Brune 345fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 346ab8d36c9SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 347ab8d36c9SPeter Brune if (flg) { 348ab8d36c9SPeter Brune ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 349fde0ff24SPeter Brune } 35087f44e3fSPeter Brune ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr); 35187f44e3fSPeter Brune if (flg) { 35287f44e3fSPeter Brune ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr); 35387f44e3fSPeter Brune } 354fde0ff24SPeter Brune 355ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 356ab8d36c9SPeter Brune if (flg) { 357ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 358ab8d36c9SPeter Brune } 359ee78dd50SPeter Brune 360928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 361928e959bSPeter 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); 362928e959bSPeter Brune if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);} 363928e959bSPeter Brune } 364928e959bSPeter Brune 36566585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 366162d76ddSPeter Brune 36766585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 368162d76ddSPeter Brune 369c8c899caSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 370c8c899caSPeter Brune if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 3710dd27c6cSPeter Brune 3720dd27c6cSPeter Brune flg = PETSC_FALSE; 3730dd27c6cSPeter Brune monflg = PETSC_TRUE; 3740dd27c6cSPeter Brune ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 3750dd27c6cSPeter Brune if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 376ab8d36c9SPeter Brune } 377ee78dd50SPeter Brune 378421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3798cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 380162d76ddSPeter Brune if (upflg) { 38166585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 382162d76ddSPeter Brune } 383162d76ddSPeter Brune if (downflg) { 38466585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 385162d76ddSPeter Brune } 386eff52c0eSPeter Brune 3879e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3889e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3899e764e56SPeter Brune if (!snes->linesearch) { 3907601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 3911a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 3929e764e56SPeter Brune } 3939e764e56SPeter Brune } 3949e764e56SPeter Brune 395ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 396ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 397ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 398421d9b32SPeter Brune PetscFunctionReturn(0); 399421d9b32SPeter Brune } 400421d9b32SPeter Brune 4019804daf3SBarry Smith #include <petscdraw.h> 402421d9b32SPeter Brune #undef __FUNCT__ 403421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 404421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 405421d9b32SPeter Brune { 406421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 407656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 408ab8d36c9SPeter Brune PetscInt i; 409421d9b32SPeter Brune PetscErrorCode ierr; 410ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 411421d9b32SPeter Brune 412421d9b32SPeter Brune PetscFunctionBegin; 413ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 414ab8d36c9SPeter Brune if (isFine) { 415251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 416656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 417421d9b32SPeter Brune if (iascii) { 418ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 419ab8d36c9SPeter Brune if (fas->galerkin) { 420ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 421421d9b32SPeter Brune } else { 422ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 423421d9b32SPeter Brune } 424ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 425ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 426ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 427ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 428ab8d36c9SPeter Brune if (!i) { 429ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 430421d9b32SPeter Brune } else { 431ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 432421d9b32SPeter Brune } 433ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 434166b3ea4SJed Brown if (smoothd) { 435ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 436166b3ea4SJed Brown } else { 437166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 438166b3ea4SJed Brown } 439ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 440ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 441ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 442ab8d36c9SPeter Brune } else if (i) { 443ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 444ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 445166b3ea4SJed Brown if (smoothu) { 446ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 447166b3ea4SJed Brown } else { 448166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 449166b3ea4SJed Brown } 450ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 451ab8d36c9SPeter Brune } 452ab8d36c9SPeter Brune } 453656ede7eSPeter Brune } else if (isdraw) { 454656ede7eSPeter Brune PetscDraw draw; 455b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 456656ede7eSPeter Brune SNES_FAS *curfas = fas; 457656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 458656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 459656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 460656ede7eSPeter Brune bottom = y - th; 461656ede7eSPeter Brune while (curfas) { 462b4375e8dSPeter Brune if (!curfas->smoothu) { 463656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 464656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 465656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 466b4375e8dSPeter Brune } else { 467b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 468b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 469b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 470b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 471b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 472b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 473b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 474b4375e8dSPeter Brune } 475656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 476656ede7eSPeter Brune bottom -= 5*th; 4771aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 4780298fd71SBarry Smith else curfas = NULL; 479656ede7eSPeter Brune } 480421d9b32SPeter Brune } 481ab8d36c9SPeter Brune } 482421d9b32SPeter Brune PetscFunctionReturn(0); 483421d9b32SPeter Brune } 484421d9b32SPeter Brune 485421d9b32SPeter Brune #undef __FUNCT__ 48691f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private" 48739bd7f45SPeter Brune /* 48839bd7f45SPeter Brune Defines the action of the downsmoother 48939bd7f45SPeter Brune */ 49091f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 491b9c2fdf1SPeter Brune { 49239bd7f45SPeter Brune PetscErrorCode ierr = 0; 493742fe5e2SPeter Brune SNESConvergedReason reason; 494ab8d36c9SPeter Brune Vec FPC; 495ab8d36c9SPeter Brune SNES smoothd; 4960dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 4976e111a19SKarl Rupp 498421d9b32SPeter Brune PetscFunctionBegin; 499ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 500e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 5010dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 502ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 5030dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 504742fe5e2SPeter Brune /* check convergence reason for the smoother */ 505ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 5061ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 507742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 508742fe5e2SPeter Brune PetscFunctionReturn(0); 509742fe5e2SPeter Brune } 5100298fd71SBarry Smith ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 5114b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 51271dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 51339bd7f45SPeter Brune PetscFunctionReturn(0); 51439bd7f45SPeter Brune } 51539bd7f45SPeter Brune 51639bd7f45SPeter Brune 51739bd7f45SPeter Brune #undef __FUNCT__ 51891f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private" 51939bd7f45SPeter Brune /* 52007144faaSPeter Brune Defines the action of the upsmoother 52139bd7f45SPeter Brune */ 5220adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 5230adebc6cSBarry Smith { 52439bd7f45SPeter Brune PetscErrorCode ierr = 0; 52539bd7f45SPeter Brune SNESConvergedReason reason; 526ab8d36c9SPeter Brune Vec FPC; 527ab8d36c9SPeter Brune SNES smoothu; 5280dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 529ab8d36c9SPeter Brune 5306e111a19SKarl Rupp PetscFunctionBegin; 531ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 5320dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 533ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 5340dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 53539bd7f45SPeter Brune /* check convergence reason for the smoother */ 536ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 5371ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 53839bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 53939bd7f45SPeter Brune PetscFunctionReturn(0); 54039bd7f45SPeter Brune } 5410298fd71SBarry Smith ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 5424b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 54371dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 54439bd7f45SPeter Brune PetscFunctionReturn(0); 54539bd7f45SPeter Brune } 54639bd7f45SPeter Brune 54739bd7f45SPeter Brune #undef __FUNCT__ 548938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 549938e4a01SJed Brown /*@ 550938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 551938e4a01SJed Brown 552938e4a01SJed Brown Collective 553938e4a01SJed Brown 554938e4a01SJed Brown Input Arguments: 555938e4a01SJed Brown . snes - SNESFAS 556938e4a01SJed Brown 557938e4a01SJed Brown Output Arguments: 558938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 559938e4a01SJed Brown 560938e4a01SJed Brown Level: developer 561938e4a01SJed Brown 562938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 563938e4a01SJed Brown @*/ 564938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 565938e4a01SJed Brown { 566938e4a01SJed Brown PetscErrorCode ierr; 567938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 568938e4a01SJed Brown 569938e4a01SJed Brown PetscFunctionBegin; 5701aa26658SKarl Rupp if (fas->rscale) { 5711aa26658SKarl Rupp ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 572f5af7f23SKarl Rupp } else if (fas->inject) { 5732a7a6963SBarry Smith ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 574ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 575938e4a01SJed Brown PetscFunctionReturn(0); 576938e4a01SJed Brown } 577938e4a01SJed Brown 578e9923e8dSJed Brown #undef __FUNCT__ 579e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 580e9923e8dSJed Brown /*@ 581e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 582e9923e8dSJed Brown 583e9923e8dSJed Brown Collective 584e9923e8dSJed Brown 585e9923e8dSJed Brown Input Arguments: 586e9923e8dSJed Brown + fine - SNES from which to restrict 587e9923e8dSJed Brown - Xfine - vector to restrict 588e9923e8dSJed Brown 589e9923e8dSJed Brown Output Arguments: 590e9923e8dSJed Brown . Xcoarse - result of restriction 591e9923e8dSJed Brown 592e9923e8dSJed Brown Level: developer 593e9923e8dSJed Brown 594e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 595e9923e8dSJed Brown @*/ 596e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 597e9923e8dSJed Brown { 598e9923e8dSJed Brown PetscErrorCode ierr; 599e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 600e9923e8dSJed Brown 601e9923e8dSJed Brown PetscFunctionBegin; 602e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 603e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 604e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 605e9923e8dSJed Brown if (fas->inject) { 606e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 607e9923e8dSJed Brown } else { 608e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 609e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 610e9923e8dSJed Brown } 611e9923e8dSJed Brown PetscFunctionReturn(0); 612e9923e8dSJed Brown } 613e9923e8dSJed Brown 614e9923e8dSJed Brown #undef __FUNCT__ 6158c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection" 61639bd7f45SPeter Brune /* 61739bd7f45SPeter Brune 61839bd7f45SPeter Brune Performs the FAS coarse correction as: 61939bd7f45SPeter Brune 620b5527d98SMatthew G. Knepley fine problem: F(x) = b 621b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c 62239bd7f45SPeter Brune 623b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 62439bd7f45SPeter Brune 62539bd7f45SPeter Brune */ 6260adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 6270adebc6cSBarry Smith { 62839bd7f45SPeter Brune PetscErrorCode ierr; 62939bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 63039bd7f45SPeter Brune SNESConvergedReason reason; 631ab8d36c9SPeter Brune SNES next; 632ab8d36c9SPeter Brune Mat restrct, interpolate; 6330dd27c6cSPeter Brune SNES_FAS *fasc; 6345fd66863SKarl Rupp 63539bd7f45SPeter Brune PetscFunctionBegin; 636ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 637ab8d36c9SPeter Brune if (next) { 6380dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 6390dd27c6cSPeter Brune 640ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 641ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 642ab8d36c9SPeter Brune 643ab8d36c9SPeter Brune X_c = next->vec_sol; 644ab8d36c9SPeter Brune Xo_c = next->work[0]; 645ab8d36c9SPeter Brune F_c = next->vec_func; 646ab8d36c9SPeter Brune B_c = next->vec_rhs; 647efe1f98aSPeter Brune 6480dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 649938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 6505609cbf2SMatthew G. Knepley /* restrict the defect: R(F(x) - b) */ 651ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 6520dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 6530dd27c6cSPeter Brune 6540dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6555609cbf2SMatthew G. Knepley /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */ 656ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 6570dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6580dd27c6cSPeter Brune 6590dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 660e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 661b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 662e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 663ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 664ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 665c90fad12SPeter Brune 666c90fad12SPeter Brune /* recurse to the next level */ 667e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 668ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 669ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 670742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 671742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 672742fe5e2SPeter Brune PetscFunctionReturn(0); 673742fe5e2SPeter Brune } 674fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 675fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 6760dd27c6cSPeter Brune 6770dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 678ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 6790dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 680293a7e31SPeter Brune } 68139bd7f45SPeter Brune PetscFunctionReturn(0); 68239bd7f45SPeter Brune } 68339bd7f45SPeter Brune 68439bd7f45SPeter Brune #undef __FUNCT__ 6852cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive" 68639bd7f45SPeter Brune /* 68739bd7f45SPeter Brune 68839bd7f45SPeter Brune The additive cycle looks like: 68939bd7f45SPeter Brune 69007144faaSPeter Brune xhat = x 69107144faaSPeter Brune xhat = dS(x, b) 69207144faaSPeter Brune x = coarsecorrection(xhat, b_d) 69307144faaSPeter Brune x = x + nu*(xhat - x); 69439bd7f45SPeter Brune (optional) x = uS(x, b) 69539bd7f45SPeter Brune 69639bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 69739bd7f45SPeter Brune 69839bd7f45SPeter Brune */ 6990adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 7000adebc6cSBarry Smith { 70107144faaSPeter Brune Vec F, B, Xhat; 70222c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 70339bd7f45SPeter Brune PetscErrorCode ierr; 70407144faaSPeter Brune SNESConvergedReason reason; 70522c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 70622c1e704SPeter Brune PetscBool lssuccess; 707ab8d36c9SPeter Brune SNES next; 708ab8d36c9SPeter Brune Mat restrct, interpolate; 7090dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 7100dd27c6cSPeter Brune 71139bd7f45SPeter Brune PetscFunctionBegin; 712ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 71339bd7f45SPeter Brune F = snes->vec_func; 71439bd7f45SPeter Brune B = snes->vec_rhs; 715e7f468e7SPeter Brune Xhat = snes->work[1]; 71607144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 71707144faaSPeter Brune /* recurse first */ 718ab8d36c9SPeter Brune if (next) { 7190dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 720ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 721ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 7220dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 72307144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 7240dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 725c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 726ab8d36c9SPeter Brune X_c = next->vec_sol; 727ab8d36c9SPeter Brune Xo_c = next->work[0]; 728ab8d36c9SPeter Brune F_c = next->vec_func; 729ab8d36c9SPeter Brune B_c = next->vec_rhs; 73039bd7f45SPeter Brune 731938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 73207144faaSPeter Brune /* restrict the defect */ 733ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 73407144faaSPeter Brune 73507144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 7360dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 737ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 7380dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 739e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 740b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 741e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 74207144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 74307144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 74407144faaSPeter Brune 74507144faaSPeter Brune /* recurse */ 746e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 747ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 74807144faaSPeter Brune 74907144faaSPeter Brune /* smooth on this level */ 75091f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 75107144faaSPeter Brune 752ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 75307144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 75407144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 75507144faaSPeter Brune PetscFunctionReturn(0); 75607144faaSPeter Brune } 75707144faaSPeter Brune 75807144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 759c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 760ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 76107144faaSPeter Brune 762ddebd997SPeter Brune /* additive correction of the coarse direction*/ 763f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 764f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 7659e764e56SPeter Brune if (!lssuccess) { 7669e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 7679e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 7689e764e56SPeter Brune PetscFunctionReturn(0); 7699e764e56SPeter Brune } 7709e764e56SPeter Brune } 771b9c2fdf1SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 77207144faaSPeter Brune } else { 77391f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 77407144faaSPeter Brune } 77539bd7f45SPeter Brune PetscFunctionReturn(0); 77639bd7f45SPeter Brune } 77739bd7f45SPeter Brune 77839bd7f45SPeter Brune #undef __FUNCT__ 7792cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative" 78039bd7f45SPeter Brune /* 78139bd7f45SPeter Brune 78239bd7f45SPeter Brune Defines the FAS cycle as: 78339bd7f45SPeter Brune 784b5527d98SMatthew G. Knepley fine problem: F(x) = b 78539bd7f45SPeter Brune coarse problem: F^c(x) = b^c 78639bd7f45SPeter Brune 787b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 78839bd7f45SPeter Brune 78939bd7f45SPeter Brune correction: 79039bd7f45SPeter Brune 79139bd7f45SPeter Brune x = x + I(x^c - Rx) 79239bd7f45SPeter Brune 79339bd7f45SPeter Brune */ 7940adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7950adebc6cSBarry Smith { 79639bd7f45SPeter Brune 79739bd7f45SPeter Brune PetscErrorCode ierr; 79839bd7f45SPeter Brune Vec F,B; 79934d65b3cSPeter Brune SNES next; 80039bd7f45SPeter Brune 80139bd7f45SPeter Brune PetscFunctionBegin; 80239bd7f45SPeter Brune F = snes->vec_func; 80339bd7f45SPeter Brune B = snes->vec_rhs; 80439bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 80534d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 80691f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 80734d65b3cSPeter Brune if (next) { 8088c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 80991f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 810fe6f9142SPeter Brune } 811fa9694d7SPeter Brune PetscFunctionReturn(0); 812421d9b32SPeter Brune } 813421d9b32SPeter Brune 814421d9b32SPeter Brune #undef __FUNCT__ 8158c94862eSPeter Brune #define __FUNCT__ "SNESFASCycleSetupPhase_Full" 8168c94862eSPeter Brune PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 8178c94862eSPeter Brune { 8188c94862eSPeter Brune SNES next; 8198c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8208c94862eSPeter Brune PetscBool isFine; 8218c94862eSPeter Brune PetscErrorCode ierr; 8228c94862eSPeter Brune 8238c94862eSPeter Brune PetscFunctionBegin; 8248c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 8258c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8268c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8278c94862eSPeter Brune fas->full_stage = 0; 8288c94862eSPeter Brune if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);} 8298c94862eSPeter Brune PetscFunctionReturn(0); 8308c94862eSPeter Brune } 8318c94862eSPeter Brune 8328c94862eSPeter Brune #undef __FUNCT__ 8338c94862eSPeter Brune #define __FUNCT__ "SNESFASCycle_Full" 8348c94862eSPeter Brune PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 8358c94862eSPeter Brune { 8368c94862eSPeter Brune PetscErrorCode ierr; 8378c94862eSPeter Brune Vec F,B; 8388c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8398c94862eSPeter Brune PetscBool isFine; 8408c94862eSPeter Brune SNES next; 8418c94862eSPeter Brune 8428c94862eSPeter Brune PetscFunctionBegin; 8438c94862eSPeter Brune F = snes->vec_func; 8448c94862eSPeter Brune B = snes->vec_rhs; 8458c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8468c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8478c94862eSPeter Brune 8488c94862eSPeter Brune if (isFine) { 8498c94862eSPeter Brune ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr); 8508c94862eSPeter Brune } 8518c94862eSPeter Brune 8528c94862eSPeter Brune if (fas->full_stage == 0) { 853928e959bSPeter Brune /* downsweep */ 8548c94862eSPeter Brune if (next) { 8558c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 856928e959bSPeter Brune if (fas->full_downsweep||isFine) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8578c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8588c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8598c94862eSPeter Brune if (fas->level != 1) next->max_its -= 1; 8608c94862eSPeter Brune } else { 8618c94862eSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8628c94862eSPeter Brune } 8638c94862eSPeter Brune fas->full_stage = 1; 8648c94862eSPeter Brune } else if (fas->full_stage == 1) { 8658c94862eSPeter Brune if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8668c94862eSPeter Brune if (next) { 8678c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8688c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8698c94862eSPeter Brune } 8708c94862eSPeter Brune } 8718c94862eSPeter Brune /* final v-cycle */ 8728c94862eSPeter Brune if (isFine) { 8738c94862eSPeter Brune if (next) { 8748c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8758c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8768c94862eSPeter Brune } 8778c94862eSPeter Brune } 8788c94862eSPeter Brune PetscFunctionReturn(0); 8798c94862eSPeter Brune } 8808c94862eSPeter Brune 8818c94862eSPeter Brune #undef __FUNCT__ 88234d65b3cSPeter Brune #define __FUNCT__ "SNESFASCycle_Kaskade" 88334d65b3cSPeter Brune PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 88434d65b3cSPeter Brune { 88534d65b3cSPeter Brune PetscErrorCode ierr; 88634d65b3cSPeter Brune Vec F,B; 88734d65b3cSPeter Brune SNES next; 88834d65b3cSPeter Brune 88934d65b3cSPeter Brune PetscFunctionBegin; 89034d65b3cSPeter Brune F = snes->vec_func; 89134d65b3cSPeter Brune B = snes->vec_rhs; 89234d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 89334d65b3cSPeter Brune if (next) { 89434d65b3cSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 89534d65b3cSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 89634d65b3cSPeter Brune } else { 89734d65b3cSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 89834d65b3cSPeter Brune } 89934d65b3cSPeter Brune PetscFunctionReturn(0); 90034d65b3cSPeter Brune } 90134d65b3cSPeter Brune 902fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE; 903fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n" 904fffbeea8SBarry Smith " title = {Composing Scalable Nonlinear Algebraic Solvers},\n" 905fffbeea8SBarry Smith " author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n" 906fffbeea8SBarry Smith " year = 2013,\n" 907fffbeea8SBarry Smith " type = Preprint,\n" 908fffbeea8SBarry Smith " number = {ANL/MCS-P2010-0112},\n" 909fffbeea8SBarry Smith " institution = {Argonne National Laboratory}\n}\n"; 910fffbeea8SBarry Smith 91134d65b3cSPeter Brune #undef __FUNCT__ 912421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 913421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 914421d9b32SPeter Brune { 915fa9694d7SPeter Brune PetscErrorCode ierr; 916fe6f9142SPeter Brune PetscInt i, maxits; 917ddb5aff1SPeter Brune Vec X, F; 918fe6f9142SPeter Brune PetscReal fnorm; 919b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 920b17ce1afSJed Brown DM dm; 921e70c42e5SPeter Brune PetscBool isFine; 922b17ce1afSJed Brown 923421d9b32SPeter Brune PetscFunctionBegin; 924fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 925fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 926fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 927fa9694d7SPeter Brune X = snes->vec_sol; 928f5a6d4f9SBarry Smith F = snes->vec_func; 929293a7e31SPeter Brune 93018a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 931293a7e31SPeter Brune /*norm setup */ 932e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 933fe6f9142SPeter Brune snes->iter = 0; 934fe6f9142SPeter Brune snes->norm = 0.; 935e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 936e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 9370dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 938fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 9390dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 940fe6f9142SPeter Brune if (snes->domainerror) { 941fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 942fe6f9142SPeter Brune PetscFunctionReturn(0); 943fe6f9142SPeter Brune } 9441aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 945e4ed7901SPeter Brune 946fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 947189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 948189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 949189a9710SBarry Smith PetscFunctionReturn(0); 950189a9710SBarry Smith } 951e4ed7901SPeter Brune 952e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 953fe6f9142SPeter Brune snes->norm = fnorm; 954e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 955a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 956fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 957fe6f9142SPeter Brune 958fe6f9142SPeter Brune /* test convergence */ 959fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 960fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 961e4ed7901SPeter Brune 962b17ce1afSJed Brown 963b9c2fdf1SPeter Brune if (isFine) { 964b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 965b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 966b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 967b17ce1afSJed Brown DM dmcoarse; 968b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 969b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 970b17ce1afSJed Brown dm = dmcoarse; 971b17ce1afSJed Brown } 972b9c2fdf1SPeter Brune } 973b17ce1afSJed Brown 974fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 975fe6f9142SPeter Brune /* Call general purpose update function */ 976646217ecSPeter Brune 977fe6f9142SPeter Brune if (snes->ops->update) { 978fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 979fe6f9142SPeter Brune } 98007144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 98191f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 9828c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_ADDITIVE) { 98391f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 9848c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_FULL) { 9858c94862eSPeter Brune ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr); 98634d65b3cSPeter Brune } else if (fas->fastype ==SNES_FAS_KASKADE) { 98734d65b3cSPeter Brune ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr); 9888c94862eSPeter Brune } else { 9894e619c20SJed Brown SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type"); 99007144faaSPeter Brune } 991742fe5e2SPeter Brune 992742fe5e2SPeter Brune /* check for FAS cycle divergence */ 9931aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 994b9c2fdf1SPeter Brune 995c90fad12SPeter Brune /* Monitor convergence */ 996e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 997c90fad12SPeter Brune snes->iter = i+1; 998e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 999a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 1000c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 1001c90fad12SPeter Brune /* Test for convergence */ 100266585501SPeter Brune if (isFine) { 1003b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1004c90fad12SPeter Brune if (snes->reason) break; 1005fe6f9142SPeter Brune } 100666585501SPeter Brune } 1007fe6f9142SPeter Brune if (i == maxits) { 1008fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1009fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1010fe6f9142SPeter Brune } 1011421d9b32SPeter Brune PetscFunctionReturn(0); 1012421d9b32SPeter Brune } 1013