1421d9b32SPeter Brune /* Defines the basic SNES object */ 2a055c8caSBarry Smith #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnes.h" I*/ 3421d9b32SPeter Brune 434d65b3cSPeter Brune const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","SNESFASType","SNES_FAS",0}; 507144faaSPeter Brune 6421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes); 7421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes); 84416b707SBarry Smith extern PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes); 9421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer); 10421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes); 11421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes); 126273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void*); 13ab8d36c9SPeter Brune extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES*); 14421d9b32SPeter Brune 151fbfccc6SJed Brown /*MC 161fbfccc6SJed Brown 171fbfccc6SJed Brown SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 181fbfccc6SJed Brown 19d3bc2379SPeter Brune The nonlinear problem is solved by correction using coarse versions 20d3bc2379SPeter Brune of the nonlinear problem. This problem is perturbed so that a projected 21d3bc2379SPeter Brune solution of the fine problem elicits no correction from the coarse problem. 22d3bc2379SPeter Brune 23d3bc2379SPeter Brune Options Database: 24d3bc2379SPeter Brune + -snes_fas_levels - The number of levels 25d3bc2379SPeter Brune . -snes_fas_cycles<1> - The number of cycles -- 1 for V, 2 for W 26928e959bSPeter Brune . -snes_fas_type<additive,multiplicative,full,kaskade> - Additive or multiplicative cycle 27d3bc2379SPeter Brune . -snes_fas_galerkin<PETSC_FALSE> - Form coarse problems by projection back upon the fine problem 28d3bc2379SPeter Brune . -snes_fas_smoothup<1> - The number of iterations of the post-smoother 29d3bc2379SPeter Brune . -snes_fas_smoothdown<1> - The number of iterations of the pre-smoother 30d3bc2379SPeter Brune . -snes_fas_monitor - Monitor progress of all of the levels 31d39b9438SMatthew G. Knepley . -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS 32d3bc2379SPeter Brune . -fas_levels_snes_ - SNES options for all smoothers 337d84e935SPeter Brune . -fas_levels_cycle_snes_ - SNES options for all cycles 34d3bc2379SPeter Brune . -fas_levels_i_snes_ - SNES options for the smoothers on level i 357d84e935SPeter Brune . -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i 36d3bc2379SPeter Brune - -fas_coarse_snes_ - SNES options for the coarsest smoother 37d3bc2379SPeter Brune 38d3bc2379SPeter Brune Notes: 39d3bc2379SPeter Brune The organization of the FAS solver is slightly different from the organization of PCMG 40d3bc2379SPeter Brune As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance. 41d3bc2379SPeter Brune The cycle SNES instance may be used for monitoring convergence on a particular level. 421fbfccc6SJed Brown 437d84e935SPeter Brune Level: beginner 441fbfccc6SJed Brown 45*4f02bc6aSBarry Smith References: 46*4f02bc6aSBarry Smith 47*4f02bc6aSBarry Smith Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 48*4f02bc6aSBarry Smith SIAM Review, 57(4), 2015 49*4f02bc6aSBarry Smith 50d3bc2379SPeter Brune .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 511fbfccc6SJed Brown M*/ 52421d9b32SPeter Brune 53421d9b32SPeter Brune #undef __FUNCT__ 54421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 558cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes) 56421d9b32SPeter Brune { 57421d9b32SPeter Brune SNES_FAS *fas; 58421d9b32SPeter Brune PetscErrorCode ierr; 59421d9b32SPeter Brune 60421d9b32SPeter Brune PetscFunctionBegin; 61421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 62421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 63421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 64421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 65421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 66421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 67421d9b32SPeter Brune 68ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 69ed020824SBarry Smith snes->usespc = PETSC_FALSE; 70ed020824SBarry Smith 7188976e71SPeter Brune if (!snes->tolerancesset) { 720e444f03SPeter Brune snes->max_funcs = 30000; 730e444f03SPeter Brune snes->max_its = 10000; 7488976e71SPeter Brune } 750e444f03SPeter Brune 76b00a9115SJed Brown ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr); 771aa26658SKarl Rupp 78421d9b32SPeter Brune snes->data = (void*) fas; 79421d9b32SPeter Brune fas->level = 0; 80293a7e31SPeter Brune fas->levels = 1; 81ee78dd50SPeter Brune fas->n_cycles = 1; 82ee78dd50SPeter Brune fas->max_up_it = 1; 83ee78dd50SPeter Brune fas->max_down_it = 1; 840298fd71SBarry Smith fas->smoothu = NULL; 850298fd71SBarry Smith fas->smoothd = NULL; 860298fd71SBarry Smith fas->next = NULL; 870298fd71SBarry Smith fas->previous = NULL; 88ab8d36c9SPeter Brune fas->fine = snes; 890298fd71SBarry Smith fas->interpolate = NULL; 900298fd71SBarry Smith fas->restrct = NULL; 910298fd71SBarry Smith fas->inject = NULL; 920298fd71SBarry Smith fas->monitor = NULL; 93cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 94ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 95928e959bSPeter Brune fas->full_downsweep = PETSC_FALSE; 960dd27c6cSPeter Brune 970298fd71SBarry Smith fas->eventsmoothsetup = 0; 980298fd71SBarry Smith fas->eventsmoothsolve = 0; 990298fd71SBarry Smith fas->eventresidual = 0; 1000298fd71SBarry Smith fas->eventinterprestrict = 0; 101efe1f98aSPeter Brune PetscFunctionReturn(0); 102efe1f98aSPeter Brune } 103efe1f98aSPeter Brune 104421d9b32SPeter Brune #undef __FUNCT__ 105421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 106421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 107421d9b32SPeter Brune { 10877df8cc4SPeter Brune PetscErrorCode ierr = 0; 109421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 110421d9b32SPeter Brune 111421d9b32SPeter Brune PetscFunctionBegin; 112ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr); 113ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 1143dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 115bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 116bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 117bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 1187cd33a7bSPeter Brune if (fas->galerkin) { 1197cd33a7bSPeter Brune ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr); 1207cd33a7bSPeter Brune ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr); 1217cd33a7bSPeter Brune } 122d477d801SPeter Brune if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);} 123421d9b32SPeter Brune PetscFunctionReturn(0); 124421d9b32SPeter Brune } 125421d9b32SPeter Brune 126421d9b32SPeter Brune #undef __FUNCT__ 127421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 128421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 129421d9b32SPeter Brune { 130421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 131742fe5e2SPeter Brune PetscErrorCode ierr = 0; 132421d9b32SPeter Brune 133421d9b32SPeter Brune PetscFunctionBegin; 134421d9b32SPeter Brune /* recursively resets and then destroys */ 13579d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 1361aa26658SKarl Rupp if (fas->next) { 1371aa26658SKarl Rupp ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 1381aa26658SKarl Rupp } 139421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 140421d9b32SPeter Brune PetscFunctionReturn(0); 141421d9b32SPeter Brune } 142421d9b32SPeter Brune 143421d9b32SPeter Brune #undef __FUNCT__ 144421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 145421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 146421d9b32SPeter Brune { 14748bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 148421d9b32SPeter Brune PetscErrorCode ierr; 149d1adcc6fSPeter Brune PetscInt dm_levels; 1503dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 151ab8d36c9SPeter Brune SNES next; 152ab8d36c9SPeter Brune PetscBool isFine; 153f89ba88eSPeter Brune SNESLineSearch linesearch; 154f89ba88eSPeter Brune SNESLineSearch slinesearch; 155f89ba88eSPeter Brune void *lsprectx,*lspostctx; 1566b2b7091SBarry Smith PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 1576b2b7091SBarry Smith PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 158eff52c0eSPeter Brune 1596b2b7091SBarry Smith PetscFunctionBegin; 160ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 161ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 162d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 163d1adcc6fSPeter Brune dm_levels++; 164cc05f883SPeter Brune if (dm_levels > fas->levels) { 1652e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 1663dccd265SPeter Brune vec_sol = snes->vec_sol; 1673dccd265SPeter Brune vec_func = snes->vec_func; 1683dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 1693dccd265SPeter Brune vec_rhs = snes->vec_rhs; 1700298fd71SBarry Smith snes->vec_sol = NULL; 1710298fd71SBarry Smith snes->vec_func = NULL; 1720298fd71SBarry Smith snes->vec_sol_update = NULL; 1730298fd71SBarry Smith snes->vec_rhs = NULL; 1743dccd265SPeter Brune 1753dccd265SPeter Brune /* reset the number of levels */ 1760298fd71SBarry Smith ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr); 177cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 1783dccd265SPeter Brune 1793dccd265SPeter Brune snes->vec_sol = vec_sol; 1803dccd265SPeter Brune snes->vec_func = vec_func; 1813dccd265SPeter Brune snes->vec_rhs = vec_rhs; 1823dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 183d1adcc6fSPeter Brune } 184d1adcc6fSPeter Brune } 185ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 186ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 1873dccd265SPeter Brune 188fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 189cc05f883SPeter Brune 190ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 191ab8d36c9SPeter Brune if (!fas->smoothd) { 192ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 193ab8d36c9SPeter Brune } 194ab8d36c9SPeter Brune 19579d9a41aSPeter Brune if (snes->dm) { 196ab8d36c9SPeter Brune /* set the smoother DMs properly */ 197ab8d36c9SPeter Brune if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 198ab8d36c9SPeter Brune ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 19979d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 200ab8d36c9SPeter Brune if (next) { 20179d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 202ab8d36c9SPeter Brune if (!next->dm) { 203ce94432eSBarry Smith ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr); 204ab8d36c9SPeter Brune ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 20579d9a41aSPeter Brune } 20679d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 20779d9a41aSPeter Brune if (!fas->interpolate) { 208ab8d36c9SPeter Brune ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 209bccf9bb3SJed Brown if (!fas->restrct) { 210bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 21179d9a41aSPeter Brune fas->restrct = fas->interpolate; 21279d9a41aSPeter Brune } 213bccf9bb3SJed Brown } 21479d9a41aSPeter Brune /* set the injection from the DM */ 21579d9a41aSPeter Brune if (!fas->inject) { 2166dbf9973SLawrence Mitchell ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr); 21779d9a41aSPeter Brune } 21879d9a41aSPeter Brune } 21979d9a41aSPeter Brune } 22079d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 22179d9a41aSPeter Brune if (fas->galerkin) { 2221aa26658SKarl Rupp if (next) { 2230298fd71SBarry Smith ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 2241aa26658SKarl Rupp } 2251aa26658SKarl Rupp if (fas->smoothd && fas->level != fas->levels - 1) { 2260298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 2271aa26658SKarl Rupp } 2281aa26658SKarl Rupp if (fas->smoothu && fas->level != fas->levels - 1) { 2290298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 2301aa26658SKarl Rupp } 23179d9a41aSPeter Brune } 23279d9a41aSPeter Brune 233534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 234534ebe21SPeter Brune if (fas->smoothd) { 235bc3f2f05SPeter Brune if (fas->level == 0 && fas->levels != 1) { 236365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); 237534ebe21SPeter Brune } else { 238365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 239534ebe21SPeter Brune } 2407fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr); 241534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); 2427601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2437601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr); 2446b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2456b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2466b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2476b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 248f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 2490dd27c6cSPeter Brune 2500dd27c6cSPeter Brune fas->smoothd->vec_sol = snes->vec_sol; 2510dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 2520dd27c6cSPeter Brune fas->smoothd->vec_sol_update = snes->vec_sol_update; 2530dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 2540dd27c6cSPeter Brune fas->smoothd->vec_func = snes->vec_func; 2550dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 2560dd27c6cSPeter Brune 2570dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2580dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr); 2590dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 260534ebe21SPeter Brune } 261534ebe21SPeter Brune 262534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 263534ebe21SPeter Brune if (fas->smoothu) { 264534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 265365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); 266534ebe21SPeter Brune } else { 267365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 268534ebe21SPeter Brune } 2697fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr); 270534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); 2717601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2727601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr); 2736b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2746b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2756b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2766b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 277f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 2780dd27c6cSPeter Brune 2790dd27c6cSPeter Brune fas->smoothu->vec_sol = snes->vec_sol; 2800dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 2810dd27c6cSPeter Brune fas->smoothu->vec_sol_update = snes->vec_sol_update; 2820dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 2830dd27c6cSPeter Brune fas->smoothu->vec_func = snes->vec_func; 2840dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 2850dd27c6cSPeter Brune 2860dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2870dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr); 2880dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 2890dd27c6cSPeter Brune 290534ebe21SPeter Brune } 291d06165b7SPeter Brune 292ab8d36c9SPeter Brune if (next) { 29379d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 294ab8d36c9SPeter Brune if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 295ab8d36c9SPeter Brune if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 2967fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); 2977601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2987601faf0SJed Brown ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr); 2996b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 3006b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 3016b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 3026b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 303f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 304ab8d36c9SPeter Brune ierr = SNESSetUp(next);CHKERRQ(ierr); 30579d9a41aSPeter Brune } 3066273346dSPeter Brune /* setup FAS work vectors */ 3076273346dSPeter Brune if (fas->galerkin) { 3086273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 3096273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 3106273346dSPeter Brune } 311421d9b32SPeter Brune PetscFunctionReturn(0); 312421d9b32SPeter Brune } 313421d9b32SPeter Brune 314421d9b32SPeter Brune #undef __FUNCT__ 315421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS" 3164416b707SBarry Smith PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes) 317421d9b32SPeter Brune { 318ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 319ee78dd50SPeter Brune PetscInt levels = 1; 32087f44e3fSPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE; 321421d9b32SPeter Brune PetscErrorCode ierr; 322ee78dd50SPeter Brune char monfilename[PETSC_MAX_PATH_LEN]; 32307144faaSPeter Brune SNESFASType fastype; 324fde0ff24SPeter Brune const char *optionsprefix; 325f1c6b773SPeter Brune SNESLineSearch linesearch; 32666585501SPeter Brune PetscInt m, n_up, n_down; 327ab8d36c9SPeter Brune SNES next; 328ab8d36c9SPeter Brune PetscBool isFine; 329421d9b32SPeter Brune 330421d9b32SPeter Brune PetscFunctionBegin; 331ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 332e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNESFAS Options-----------------------------------");CHKERRQ(ierr); 333ee78dd50SPeter Brune 334ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 335ab8d36c9SPeter Brune if (isFine) { 336ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 337c732cbdbSBarry Smith if (!flg && snes->dm) { 338c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 339c732cbdbSBarry Smith levels++; 340d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 341c732cbdbSBarry Smith } 3420298fd71SBarry Smith ierr = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr); 34307144faaSPeter Brune fastype = fas->fastype; 34407144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 34507144faaSPeter Brune if (flg) { 34607144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 34707144faaSPeter Brune } 348ee78dd50SPeter Brune 349fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 350ab8d36c9SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 351ab8d36c9SPeter Brune if (flg) { 352ab8d36c9SPeter Brune ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 353fde0ff24SPeter Brune } 35487f44e3fSPeter Brune ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr); 35587f44e3fSPeter Brune if (flg) { 35687f44e3fSPeter Brune ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr); 35787f44e3fSPeter Brune } 358fde0ff24SPeter Brune 359ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 360ab8d36c9SPeter Brune if (flg) { 361ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 362ab8d36c9SPeter Brune } 363ee78dd50SPeter Brune 364928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 365928e959bSPeter 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); 366928e959bSPeter Brune if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);} 367928e959bSPeter Brune } 368928e959bSPeter Brune 36966585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 370162d76ddSPeter Brune 37166585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 372162d76ddSPeter Brune 373c8c899caSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 374c8c899caSPeter Brune if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 3750dd27c6cSPeter Brune 3760dd27c6cSPeter Brune flg = PETSC_FALSE; 3770dd27c6cSPeter Brune monflg = PETSC_TRUE; 3780dd27c6cSPeter Brune ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 3790dd27c6cSPeter Brune if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 380ab8d36c9SPeter Brune } 381ee78dd50SPeter Brune 382421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3838cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 384162d76ddSPeter Brune if (upflg) { 38566585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 386162d76ddSPeter Brune } 387162d76ddSPeter Brune if (downflg) { 38866585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 389162d76ddSPeter Brune } 390eff52c0eSPeter Brune 3919e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3929e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3939e764e56SPeter Brune if (!snes->linesearch) { 3947601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 3951a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 3969e764e56SPeter Brune } 3979e764e56SPeter Brune } 3989e764e56SPeter Brune 399ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 400ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 401ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 402421d9b32SPeter Brune PetscFunctionReturn(0); 403421d9b32SPeter Brune } 404421d9b32SPeter Brune 4059804daf3SBarry Smith #include <petscdraw.h> 406421d9b32SPeter Brune #undef __FUNCT__ 407421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 408421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 409421d9b32SPeter Brune { 410421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 411656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 412ab8d36c9SPeter Brune PetscInt i; 413421d9b32SPeter Brune PetscErrorCode ierr; 414ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 415421d9b32SPeter Brune 416421d9b32SPeter Brune PetscFunctionBegin; 417ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 418ab8d36c9SPeter Brune if (isFine) { 419251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 420656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 421421d9b32SPeter Brune if (iascii) { 422ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 423ab8d36c9SPeter Brune if (fas->galerkin) { 424ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 425421d9b32SPeter Brune } else { 426ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 427421d9b32SPeter Brune } 428ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 429ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 430ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 431ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 432ab8d36c9SPeter Brune if (!i) { 433ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 434421d9b32SPeter Brune } else { 435ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 436421d9b32SPeter Brune } 437ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 438166b3ea4SJed Brown if (smoothd) { 439ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 440166b3ea4SJed Brown } else { 441166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 442166b3ea4SJed Brown } 443ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 444ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 445ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 446ab8d36c9SPeter Brune } else if (i) { 447ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 448ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 449166b3ea4SJed Brown if (smoothu) { 450ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 451166b3ea4SJed Brown } else { 452166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 453166b3ea4SJed Brown } 454ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 455ab8d36c9SPeter Brune } 456ab8d36c9SPeter Brune } 457656ede7eSPeter Brune } else if (isdraw) { 458656ede7eSPeter Brune PetscDraw draw; 459b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 460656ede7eSPeter Brune SNES_FAS *curfas = fas; 461656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 462656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 463656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 464656ede7eSPeter Brune bottom = y - th; 465656ede7eSPeter Brune while (curfas) { 466b4375e8dSPeter Brune if (!curfas->smoothu) { 467656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 468656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 469656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 470b4375e8dSPeter Brune } else { 471b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 472b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 473b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 474b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 475b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 476b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 477b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 478b4375e8dSPeter Brune } 479656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 480656ede7eSPeter Brune bottom -= 5*th; 4811aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 4820298fd71SBarry Smith else curfas = NULL; 483656ede7eSPeter Brune } 484421d9b32SPeter Brune } 485ab8d36c9SPeter Brune } 486421d9b32SPeter Brune PetscFunctionReturn(0); 487421d9b32SPeter Brune } 488421d9b32SPeter Brune 489421d9b32SPeter Brune #undef __FUNCT__ 49091f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private" 49139bd7f45SPeter Brune /* 49239bd7f45SPeter Brune Defines the action of the downsmoother 49339bd7f45SPeter Brune */ 49491f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 495b9c2fdf1SPeter Brune { 49639bd7f45SPeter Brune PetscErrorCode ierr = 0; 497742fe5e2SPeter Brune SNESConvergedReason reason; 498ab8d36c9SPeter Brune Vec FPC; 499ab8d36c9SPeter Brune SNES smoothd; 5000dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 5016e111a19SKarl Rupp 502421d9b32SPeter Brune PetscFunctionBegin; 503ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 504e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 5050dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 506ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 5070dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 508742fe5e2SPeter Brune /* check convergence reason for the smoother */ 509ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 5101ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 511742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 512742fe5e2SPeter Brune PetscFunctionReturn(0); 513742fe5e2SPeter Brune } 5140298fd71SBarry Smith ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 5154b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 51671dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 51739bd7f45SPeter Brune PetscFunctionReturn(0); 51839bd7f45SPeter Brune } 51939bd7f45SPeter Brune 52039bd7f45SPeter Brune 52139bd7f45SPeter Brune #undef __FUNCT__ 52291f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private" 52339bd7f45SPeter Brune /* 52407144faaSPeter Brune Defines the action of the upsmoother 52539bd7f45SPeter Brune */ 5260adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 5270adebc6cSBarry Smith { 52839bd7f45SPeter Brune PetscErrorCode ierr = 0; 52939bd7f45SPeter Brune SNESConvergedReason reason; 530ab8d36c9SPeter Brune Vec FPC; 531ab8d36c9SPeter Brune SNES smoothu; 5320dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 533ab8d36c9SPeter Brune 5346e111a19SKarl Rupp PetscFunctionBegin; 535ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 5360dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 537ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 5380dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 53939bd7f45SPeter Brune /* check convergence reason for the smoother */ 540ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 5411ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 54239bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 54339bd7f45SPeter Brune PetscFunctionReturn(0); 54439bd7f45SPeter Brune } 5450298fd71SBarry Smith ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 5464b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 54771dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 54839bd7f45SPeter Brune PetscFunctionReturn(0); 54939bd7f45SPeter Brune } 55039bd7f45SPeter Brune 55139bd7f45SPeter Brune #undef __FUNCT__ 552938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 553938e4a01SJed Brown /*@ 554938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 555938e4a01SJed Brown 556938e4a01SJed Brown Collective 557938e4a01SJed Brown 558938e4a01SJed Brown Input Arguments: 559938e4a01SJed Brown . snes - SNESFAS 560938e4a01SJed Brown 561938e4a01SJed Brown Output Arguments: 562938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 563938e4a01SJed Brown 564938e4a01SJed Brown Level: developer 565938e4a01SJed Brown 566938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 567938e4a01SJed Brown @*/ 568938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 569938e4a01SJed Brown { 570938e4a01SJed Brown PetscErrorCode ierr; 571938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 572938e4a01SJed Brown 573938e4a01SJed Brown PetscFunctionBegin; 5741aa26658SKarl Rupp if (fas->rscale) { 5751aa26658SKarl Rupp ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 576f5af7f23SKarl Rupp } else if (fas->inject) { 5772a7a6963SBarry Smith ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 578ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 579938e4a01SJed Brown PetscFunctionReturn(0); 580938e4a01SJed Brown } 581938e4a01SJed Brown 582e9923e8dSJed Brown #undef __FUNCT__ 583e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 584e9923e8dSJed Brown /*@ 585e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 586e9923e8dSJed Brown 587e9923e8dSJed Brown Collective 588e9923e8dSJed Brown 589e9923e8dSJed Brown Input Arguments: 590e9923e8dSJed Brown + fine - SNES from which to restrict 591e9923e8dSJed Brown - Xfine - vector to restrict 592e9923e8dSJed Brown 593e9923e8dSJed Brown Output Arguments: 594e9923e8dSJed Brown . Xcoarse - result of restriction 595e9923e8dSJed Brown 596e9923e8dSJed Brown Level: developer 597e9923e8dSJed Brown 598e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 599e9923e8dSJed Brown @*/ 600e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 601e9923e8dSJed Brown { 602e9923e8dSJed Brown PetscErrorCode ierr; 603e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 604e9923e8dSJed Brown 605e9923e8dSJed Brown PetscFunctionBegin; 606e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 607e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 608e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 609e9923e8dSJed Brown if (fas->inject) { 610e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 611e9923e8dSJed Brown } else { 612e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 613e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 614e9923e8dSJed Brown } 615e9923e8dSJed Brown PetscFunctionReturn(0); 616e9923e8dSJed Brown } 617e9923e8dSJed Brown 618e9923e8dSJed Brown #undef __FUNCT__ 6198c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection" 62039bd7f45SPeter Brune /* 62139bd7f45SPeter Brune 62239bd7f45SPeter Brune Performs the FAS coarse correction as: 62339bd7f45SPeter Brune 624b5527d98SMatthew G. Knepley fine problem: F(x) = b 625b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c 62639bd7f45SPeter Brune 627b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 62839bd7f45SPeter Brune 62939bd7f45SPeter Brune */ 6300adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 6310adebc6cSBarry Smith { 63239bd7f45SPeter Brune PetscErrorCode ierr; 63339bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 63439bd7f45SPeter Brune SNESConvergedReason reason; 635ab8d36c9SPeter Brune SNES next; 636ab8d36c9SPeter Brune Mat restrct, interpolate; 6370dd27c6cSPeter Brune SNES_FAS *fasc; 6385fd66863SKarl Rupp 63939bd7f45SPeter Brune PetscFunctionBegin; 640ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 641ab8d36c9SPeter Brune if (next) { 6420dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 6430dd27c6cSPeter Brune 644ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 645ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 646ab8d36c9SPeter Brune 647ab8d36c9SPeter Brune X_c = next->vec_sol; 648ab8d36c9SPeter Brune Xo_c = next->work[0]; 649ab8d36c9SPeter Brune F_c = next->vec_func; 650ab8d36c9SPeter Brune B_c = next->vec_rhs; 651efe1f98aSPeter Brune 6520dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 653938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 6545609cbf2SMatthew G. Knepley /* restrict the defect: R(F(x) - b) */ 655ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 6560dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 6570dd27c6cSPeter Brune 6580dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6595609cbf2SMatthew G. Knepley /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */ 660ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 6610dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6620dd27c6cSPeter Brune 6630dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 664e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 665b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 666e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 667ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 668ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 669c90fad12SPeter Brune 670c90fad12SPeter Brune /* recurse to the next level */ 671e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 672ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 673ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 674742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 675742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 676742fe5e2SPeter Brune PetscFunctionReturn(0); 677742fe5e2SPeter Brune } 678fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 679fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 6800dd27c6cSPeter Brune 6810dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 682ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 6830dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 684293a7e31SPeter Brune } 68539bd7f45SPeter Brune PetscFunctionReturn(0); 68639bd7f45SPeter Brune } 68739bd7f45SPeter Brune 68839bd7f45SPeter Brune #undef __FUNCT__ 6892cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive" 69039bd7f45SPeter Brune /* 69139bd7f45SPeter Brune 69239bd7f45SPeter Brune The additive cycle looks like: 69339bd7f45SPeter Brune 69407144faaSPeter Brune xhat = x 69507144faaSPeter Brune xhat = dS(x, b) 69607144faaSPeter Brune x = coarsecorrection(xhat, b_d) 69707144faaSPeter Brune x = x + nu*(xhat - x); 69839bd7f45SPeter Brune (optional) x = uS(x, b) 69939bd7f45SPeter Brune 70039bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 70139bd7f45SPeter Brune 70239bd7f45SPeter Brune */ 7030adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 7040adebc6cSBarry Smith { 70507144faaSPeter Brune Vec F, B, Xhat; 70622c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 70739bd7f45SPeter Brune PetscErrorCode ierr; 70807144faaSPeter Brune SNESConvergedReason reason; 70922c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 710422a814eSBarry Smith SNESLineSearchReason lsresult; 711ab8d36c9SPeter Brune SNES next; 712ab8d36c9SPeter Brune Mat restrct, interpolate; 7130dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 7140dd27c6cSPeter Brune 71539bd7f45SPeter Brune PetscFunctionBegin; 716ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 71739bd7f45SPeter Brune F = snes->vec_func; 71839bd7f45SPeter Brune B = snes->vec_rhs; 719e7f468e7SPeter Brune Xhat = snes->work[1]; 72007144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 72107144faaSPeter Brune /* recurse first */ 722ab8d36c9SPeter Brune if (next) { 7230dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 724ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 725ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 7260dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 72707144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 7280dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 729c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 730ab8d36c9SPeter Brune X_c = next->vec_sol; 731ab8d36c9SPeter Brune Xo_c = next->work[0]; 732ab8d36c9SPeter Brune F_c = next->vec_func; 733ab8d36c9SPeter Brune B_c = next->vec_rhs; 73439bd7f45SPeter Brune 735938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 73607144faaSPeter Brune /* restrict the defect */ 737ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 73807144faaSPeter Brune 73907144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 7400dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 741ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 7420dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 743e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 744b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 745e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 74607144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 74707144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 74807144faaSPeter Brune 74907144faaSPeter Brune /* recurse */ 750e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 751ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 75207144faaSPeter Brune 75307144faaSPeter Brune /* smooth on this level */ 75491f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 75507144faaSPeter Brune 756ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 75707144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 75807144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 75907144faaSPeter Brune PetscFunctionReturn(0); 76007144faaSPeter Brune } 76107144faaSPeter Brune 76207144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 763c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 764ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 76507144faaSPeter Brune 766ddebd997SPeter Brune /* additive correction of the coarse direction*/ 767f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 768422a814eSBarry Smith ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr); 769422a814eSBarry Smith ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 770422a814eSBarry Smith if (lsresult) { 7719e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 7729e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 7739e764e56SPeter Brune PetscFunctionReturn(0); 7749e764e56SPeter Brune } 7759e764e56SPeter Brune } 77607144faaSPeter Brune } else { 77791f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 77807144faaSPeter Brune } 77939bd7f45SPeter Brune PetscFunctionReturn(0); 78039bd7f45SPeter Brune } 78139bd7f45SPeter Brune 78239bd7f45SPeter Brune #undef __FUNCT__ 7832cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative" 78439bd7f45SPeter Brune /* 78539bd7f45SPeter Brune 78639bd7f45SPeter Brune Defines the FAS cycle as: 78739bd7f45SPeter Brune 788b5527d98SMatthew G. Knepley fine problem: F(x) = b 78939bd7f45SPeter Brune coarse problem: F^c(x) = b^c 79039bd7f45SPeter Brune 791b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 79239bd7f45SPeter Brune 79339bd7f45SPeter Brune correction: 79439bd7f45SPeter Brune 79539bd7f45SPeter Brune x = x + I(x^c - Rx) 79639bd7f45SPeter Brune 79739bd7f45SPeter Brune */ 7980adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7990adebc6cSBarry Smith { 80039bd7f45SPeter Brune 80139bd7f45SPeter Brune PetscErrorCode ierr; 80239bd7f45SPeter Brune Vec F,B; 80334d65b3cSPeter Brune SNES next; 80439bd7f45SPeter Brune 80539bd7f45SPeter Brune PetscFunctionBegin; 80639bd7f45SPeter Brune F = snes->vec_func; 80739bd7f45SPeter Brune B = snes->vec_rhs; 80839bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 80934d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 81091f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 81134d65b3cSPeter Brune if (next) { 8128c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 81391f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 814fe6f9142SPeter Brune } 815fa9694d7SPeter Brune PetscFunctionReturn(0); 816421d9b32SPeter Brune } 817421d9b32SPeter Brune 818421d9b32SPeter Brune #undef __FUNCT__ 8198c94862eSPeter Brune #define __FUNCT__ "SNESFASCycleSetupPhase_Full" 8208c94862eSPeter Brune PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 8218c94862eSPeter Brune { 8228c94862eSPeter Brune SNES next; 8238c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8248c94862eSPeter Brune PetscBool isFine; 8258c94862eSPeter Brune PetscErrorCode ierr; 8268c94862eSPeter Brune 8278c94862eSPeter Brune PetscFunctionBegin; 8288c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 8298c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8308c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8318c94862eSPeter Brune fas->full_stage = 0; 8328c94862eSPeter Brune if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);} 8338c94862eSPeter Brune PetscFunctionReturn(0); 8348c94862eSPeter Brune } 8358c94862eSPeter Brune 8368c94862eSPeter Brune #undef __FUNCT__ 8378c94862eSPeter Brune #define __FUNCT__ "SNESFASCycle_Full" 8388c94862eSPeter Brune PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 8398c94862eSPeter Brune { 8408c94862eSPeter Brune PetscErrorCode ierr; 8418c94862eSPeter Brune Vec F,B; 8428c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 8438c94862eSPeter Brune PetscBool isFine; 8448c94862eSPeter Brune SNES next; 8458c94862eSPeter Brune 8468c94862eSPeter Brune PetscFunctionBegin; 8478c94862eSPeter Brune F = snes->vec_func; 8488c94862eSPeter Brune B = snes->vec_rhs; 8498c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8508c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8518c94862eSPeter Brune 8528c94862eSPeter Brune if (isFine) { 8538c94862eSPeter Brune ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr); 8548c94862eSPeter Brune } 8558c94862eSPeter Brune 8568c94862eSPeter Brune if (fas->full_stage == 0) { 857928e959bSPeter Brune /* downsweep */ 8588c94862eSPeter Brune if (next) { 8598c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 860928e959bSPeter Brune if (fas->full_downsweep||isFine) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8618c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8628c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8638c94862eSPeter Brune if (fas->level != 1) next->max_its -= 1; 8648c94862eSPeter Brune } else { 865a3a80b83SMatthew G. Knepley /* The smoother on the coarse level is the coarse solver */ 8668c94862eSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8678c94862eSPeter Brune } 8688c94862eSPeter Brune fas->full_stage = 1; 8698c94862eSPeter Brune } else if (fas->full_stage == 1) { 8708c94862eSPeter Brune if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8718c94862eSPeter Brune if (next) { 8728c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8738c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8748c94862eSPeter Brune } 8758c94862eSPeter Brune } 8768c94862eSPeter Brune /* final v-cycle */ 8778c94862eSPeter Brune if (isFine) { 8788c94862eSPeter Brune if (next) { 8798c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8808c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8818c94862eSPeter Brune } 8828c94862eSPeter Brune } 8838c94862eSPeter Brune PetscFunctionReturn(0); 8848c94862eSPeter Brune } 8858c94862eSPeter Brune 8868c94862eSPeter Brune #undef __FUNCT__ 88734d65b3cSPeter Brune #define __FUNCT__ "SNESFASCycle_Kaskade" 88834d65b3cSPeter Brune PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 88934d65b3cSPeter Brune { 89034d65b3cSPeter Brune PetscErrorCode ierr; 89134d65b3cSPeter Brune Vec F,B; 89234d65b3cSPeter Brune SNES next; 89334d65b3cSPeter Brune 89434d65b3cSPeter Brune PetscFunctionBegin; 89534d65b3cSPeter Brune F = snes->vec_func; 89634d65b3cSPeter Brune B = snes->vec_rhs; 89734d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 89834d65b3cSPeter Brune if (next) { 89934d65b3cSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 90034d65b3cSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 90134d65b3cSPeter Brune } else { 90234d65b3cSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 90334d65b3cSPeter Brune } 90434d65b3cSPeter Brune PetscFunctionReturn(0); 90534d65b3cSPeter Brune } 90634d65b3cSPeter Brune 907fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE; 908fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n" 909fffbeea8SBarry Smith " title = {Composing Scalable Nonlinear Algebraic Solvers},\n" 910fffbeea8SBarry Smith " author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n" 911fffbeea8SBarry Smith " year = 2013,\n" 912fffbeea8SBarry Smith " type = Preprint,\n" 913fffbeea8SBarry Smith " number = {ANL/MCS-P2010-0112},\n" 914fffbeea8SBarry Smith " institution = {Argonne National Laboratory}\n}\n"; 915fffbeea8SBarry Smith 91634d65b3cSPeter Brune #undef __FUNCT__ 917421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 918421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 919421d9b32SPeter Brune { 920fa9694d7SPeter Brune PetscErrorCode ierr; 921fe6f9142SPeter Brune PetscInt i, maxits; 922ddb5aff1SPeter Brune Vec X, F; 923fe6f9142SPeter Brune PetscReal fnorm; 924b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 925b17ce1afSJed Brown DM dm; 926e70c42e5SPeter Brune PetscBool isFine; 927b17ce1afSJed Brown 928421d9b32SPeter Brune PetscFunctionBegin; 929c579b300SPatrick Farrell 9306c4ed002SBarry Smith if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 931c579b300SPatrick Farrell 932fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 933fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 934fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 935fa9694d7SPeter Brune X = snes->vec_sol; 936f5a6d4f9SBarry Smith F = snes->vec_func; 937293a7e31SPeter Brune 93818a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 939293a7e31SPeter Brune /*norm setup */ 940e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 941fe6f9142SPeter Brune snes->iter = 0; 942fe6f9142SPeter Brune snes->norm = 0.; 943e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 944e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 9450dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 946fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 9470dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 9481aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 949e4ed7901SPeter Brune 950fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 951422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 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); 9886c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type"); 989742fe5e2SPeter Brune 990742fe5e2SPeter Brune /* check for FAS cycle divergence */ 9911aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 992b9c2fdf1SPeter Brune 993c90fad12SPeter Brune /* Monitor convergence */ 994e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 995c90fad12SPeter Brune snes->iter = i+1; 996e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 997a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 998c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 999c90fad12SPeter Brune /* Test for convergence */ 100066585501SPeter Brune if (isFine) { 1001b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 1002c90fad12SPeter Brune if (snes->reason) break; 1003fe6f9142SPeter Brune } 100466585501SPeter Brune } 1005fe6f9142SPeter Brune if (i == maxits) { 1006fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 1007fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 1008fe6f9142SPeter Brune } 1009421d9b32SPeter Brune PetscFunctionReturn(0); 1010421d9b32SPeter Brune } 1011