1421d9b32SPeter Brune /* Defines the basic SNES object */ 26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnesfas.h" I*/ 3421d9b32SPeter Brune 46a6fc655SJed Brown const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","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 26d3bc2379SPeter Brune . -snes_fas_type<additive, multiplicative> - 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 31d3bc2379SPeter Brune . -fas_levels_snes_ - SNES options for all smoothers 327d84e935SPeter Brune . -fas_levels_cycle_snes_ - SNES options for all cycles 33d3bc2379SPeter Brune . -fas_levels_i_snes_ - SNES options for the smoothers on level i 347d84e935SPeter Brune . -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i 35d3bc2379SPeter Brune - -fas_coarse_snes_ - SNES options for the coarsest smoother 36d3bc2379SPeter Brune 37d3bc2379SPeter Brune Notes: 38d3bc2379SPeter Brune The organization of the FAS solver is slightly different from the organization of PCMG 39d3bc2379SPeter Brune As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance. 40d3bc2379SPeter Brune The cycle SNES instance may be used for monitoring convergence on a particular level. 411fbfccc6SJed Brown 427d84e935SPeter Brune Level: beginner 431fbfccc6SJed Brown 44d3bc2379SPeter Brune .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 451fbfccc6SJed Brown M*/ 46421d9b32SPeter Brune 47421d9b32SPeter Brune #undef __FUNCT__ 48421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS" 498cc058d9SJed Brown PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes) 50421d9b32SPeter Brune { 51421d9b32SPeter Brune SNES_FAS *fas; 52421d9b32SPeter Brune PetscErrorCode ierr; 53421d9b32SPeter Brune 54421d9b32SPeter Brune PetscFunctionBegin; 55421d9b32SPeter Brune snes->ops->destroy = SNESDestroy_FAS; 56421d9b32SPeter Brune snes->ops->setup = SNESSetUp_FAS; 57421d9b32SPeter Brune snes->ops->setfromoptions = SNESSetFromOptions_FAS; 58421d9b32SPeter Brune snes->ops->view = SNESView_FAS; 59421d9b32SPeter Brune snes->ops->solve = SNESSolve_FAS; 60421d9b32SPeter Brune snes->ops->reset = SNESReset_FAS; 61421d9b32SPeter Brune 62ed020824SBarry Smith snes->usesksp = PETSC_FALSE; 63ed020824SBarry Smith snes->usespc = PETSC_FALSE; 64ed020824SBarry Smith 6588976e71SPeter Brune if (!snes->tolerancesset) { 660e444f03SPeter Brune snes->max_funcs = 30000; 670e444f03SPeter Brune snes->max_its = 10000; 6888976e71SPeter Brune } 690e444f03SPeter Brune 70421d9b32SPeter Brune ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr); 711aa26658SKarl Rupp 72421d9b32SPeter Brune snes->data = (void*) fas; 73421d9b32SPeter Brune fas->level = 0; 74293a7e31SPeter Brune fas->levels = 1; 75ee78dd50SPeter Brune fas->n_cycles = 1; 76ee78dd50SPeter Brune fas->max_up_it = 1; 77ee78dd50SPeter Brune fas->max_down_it = 1; 780298fd71SBarry Smith fas->smoothu = NULL; 790298fd71SBarry Smith fas->smoothd = NULL; 800298fd71SBarry Smith fas->next = NULL; 810298fd71SBarry Smith fas->previous = NULL; 82ab8d36c9SPeter Brune fas->fine = snes; 830298fd71SBarry Smith fas->interpolate = NULL; 840298fd71SBarry Smith fas->restrct = NULL; 850298fd71SBarry Smith fas->inject = NULL; 860298fd71SBarry Smith fas->monitor = NULL; 87cc05f883SPeter Brune fas->usedmfornumberoflevels = PETSC_FALSE; 88ddebd997SPeter Brune fas->fastype = SNES_FAS_MULTIPLICATIVE; 890dd27c6cSPeter Brune 900298fd71SBarry Smith fas->eventsmoothsetup = 0; 910298fd71SBarry Smith fas->eventsmoothsolve = 0; 920298fd71SBarry Smith fas->eventresidual = 0; 930298fd71SBarry Smith fas->eventinterprestrict = 0; 94efe1f98aSPeter Brune PetscFunctionReturn(0); 95efe1f98aSPeter Brune } 96efe1f98aSPeter Brune 97421d9b32SPeter Brune #undef __FUNCT__ 98421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS" 99421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes) 100421d9b32SPeter Brune { 10177df8cc4SPeter Brune PetscErrorCode ierr = 0; 102421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 103421d9b32SPeter Brune 104421d9b32SPeter Brune PetscFunctionBegin; 105ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr); 106ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 1073dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 108bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 109bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 110bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 1117cd33a7bSPeter Brune if (fas->galerkin) { 1127cd33a7bSPeter Brune ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr); 1137cd33a7bSPeter Brune ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr); 1147cd33a7bSPeter Brune } 115*d477d801SPeter Brune if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);} 116421d9b32SPeter Brune PetscFunctionReturn(0); 117421d9b32SPeter Brune } 118421d9b32SPeter Brune 119421d9b32SPeter Brune #undef __FUNCT__ 120421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS" 121421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes) 122421d9b32SPeter Brune { 123421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 124742fe5e2SPeter Brune PetscErrorCode ierr = 0; 125421d9b32SPeter Brune 126421d9b32SPeter Brune PetscFunctionBegin; 127421d9b32SPeter Brune /* recursively resets and then destroys */ 12879d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 1291aa26658SKarl Rupp if (fas->next) { 1301aa26658SKarl Rupp ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 1311aa26658SKarl Rupp } 132421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 133421d9b32SPeter Brune PetscFunctionReturn(0); 134421d9b32SPeter Brune } 135421d9b32SPeter Brune 136421d9b32SPeter Brune #undef __FUNCT__ 137421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS" 138421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes) 139421d9b32SPeter Brune { 14048bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 141421d9b32SPeter Brune PetscErrorCode ierr; 142efe1f98aSPeter Brune VecScatter injscatter; 143d1adcc6fSPeter Brune PetscInt dm_levels; 1443dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 145ab8d36c9SPeter Brune SNES next; 146ab8d36c9SPeter Brune PetscBool isFine; 147f89ba88eSPeter Brune SNESLineSearch linesearch; 148f89ba88eSPeter Brune SNESLineSearch slinesearch; 149f89ba88eSPeter Brune void *lsprectx,*lspostctx; 1506b2b7091SBarry Smith PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 1516b2b7091SBarry Smith PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 152eff52c0eSPeter Brune 1536b2b7091SBarry Smith PetscFunctionBegin; 154ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 155ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 156d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 157d1adcc6fSPeter Brune dm_levels++; 158cc05f883SPeter Brune if (dm_levels > fas->levels) { 1592e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 1603dccd265SPeter Brune vec_sol = snes->vec_sol; 1613dccd265SPeter Brune vec_func = snes->vec_func; 1623dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 1633dccd265SPeter Brune vec_rhs = snes->vec_rhs; 1640298fd71SBarry Smith snes->vec_sol = NULL; 1650298fd71SBarry Smith snes->vec_func = NULL; 1660298fd71SBarry Smith snes->vec_sol_update = NULL; 1670298fd71SBarry Smith snes->vec_rhs = NULL; 1683dccd265SPeter Brune 1693dccd265SPeter Brune /* reset the number of levels */ 1700298fd71SBarry Smith ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr); 171cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 1723dccd265SPeter Brune 1733dccd265SPeter Brune snes->vec_sol = vec_sol; 1743dccd265SPeter Brune snes->vec_func = vec_func; 1753dccd265SPeter Brune snes->vec_rhs = vec_rhs; 1763dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 177d1adcc6fSPeter Brune } 178d1adcc6fSPeter Brune } 179ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 180ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 1813dccd265SPeter Brune 182fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 183cc05f883SPeter Brune 184ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 185ab8d36c9SPeter Brune if (!fas->smoothd) { 186ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 187ab8d36c9SPeter Brune } 188ab8d36c9SPeter Brune 18979d9a41aSPeter Brune if (snes->dm) { 190ab8d36c9SPeter Brune /* set the smoother DMs properly */ 191ab8d36c9SPeter Brune if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 192ab8d36c9SPeter Brune ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 19379d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 194ab8d36c9SPeter Brune if (next) { 19579d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 196ab8d36c9SPeter Brune if (!next->dm) { 197ce94432eSBarry Smith ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr); 198ab8d36c9SPeter Brune ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 19979d9a41aSPeter Brune } 20079d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 20179d9a41aSPeter Brune if (!fas->interpolate) { 202ab8d36c9SPeter Brune ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 203bccf9bb3SJed Brown if (!fas->restrct) { 204bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 20579d9a41aSPeter Brune fas->restrct = fas->interpolate; 20679d9a41aSPeter Brune } 207bccf9bb3SJed Brown } 20879d9a41aSPeter Brune /* set the injection from the DM */ 20979d9a41aSPeter Brune if (!fas->inject) { 210ab8d36c9SPeter Brune ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr); 211ce94432eSBarry Smith ierr = MatCreateScatter(PetscObjectComm((PetscObject)snes), injscatter, &fas->inject);CHKERRQ(ierr); 21279d9a41aSPeter Brune ierr = VecScatterDestroy(&injscatter);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; 3164d26bfa5SPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = 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 } 350fde0ff24SPeter Brune 351ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 352ab8d36c9SPeter Brune if (flg) { 353ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 354ab8d36c9SPeter Brune } 355ee78dd50SPeter Brune 35666585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 357162d76ddSPeter Brune 35866585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 359162d76ddSPeter Brune 360c8c899caSPeter Brune ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr); 361c8c899caSPeter Brune if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr); 3620dd27c6cSPeter Brune 3630dd27c6cSPeter Brune flg = PETSC_FALSE; 3640dd27c6cSPeter Brune monflg = PETSC_TRUE; 3650dd27c6cSPeter Brune ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 3660dd27c6cSPeter Brune if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 367ab8d36c9SPeter Brune } 368ee78dd50SPeter Brune 369421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 3708cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 371162d76ddSPeter Brune if (upflg) { 37266585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 373162d76ddSPeter Brune } 374162d76ddSPeter Brune if (downflg) { 37566585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 376162d76ddSPeter Brune } 377eff52c0eSPeter Brune 3789e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3799e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3809e764e56SPeter Brune if (!snes->linesearch) { 3817601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 3821a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 3839e764e56SPeter Brune } 3849e764e56SPeter Brune } 3859e764e56SPeter Brune 386ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 387ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 388ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 389421d9b32SPeter Brune PetscFunctionReturn(0); 390421d9b32SPeter Brune } 391421d9b32SPeter Brune 3929804daf3SBarry Smith #include <petscdraw.h> 393421d9b32SPeter Brune #undef __FUNCT__ 394421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS" 395421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 396421d9b32SPeter Brune { 397421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 398656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 399ab8d36c9SPeter Brune PetscInt i; 400421d9b32SPeter Brune PetscErrorCode ierr; 401ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 402421d9b32SPeter Brune 403421d9b32SPeter Brune PetscFunctionBegin; 404ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 405ab8d36c9SPeter Brune if (isFine) { 406251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 407656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 408421d9b32SPeter Brune if (iascii) { 409ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 410ab8d36c9SPeter Brune if (fas->galerkin) { 411ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 412421d9b32SPeter Brune } else { 413ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 414421d9b32SPeter Brune } 415ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 416ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 417ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 418ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 419ab8d36c9SPeter Brune if (!i) { 420ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 421421d9b32SPeter Brune } else { 422ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 423421d9b32SPeter Brune } 424ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 425166b3ea4SJed Brown if (smoothd) { 426ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 427166b3ea4SJed Brown } else { 428166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 429166b3ea4SJed Brown } 430ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 431ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 432ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 433ab8d36c9SPeter Brune } else if (i) { 434ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 435ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 436166b3ea4SJed Brown if (smoothu) { 437ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 438166b3ea4SJed Brown } else { 439166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 440166b3ea4SJed Brown } 441ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 442ab8d36c9SPeter Brune } 443ab8d36c9SPeter Brune } 444656ede7eSPeter Brune } else if (isdraw) { 445656ede7eSPeter Brune PetscDraw draw; 446b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 447656ede7eSPeter Brune SNES_FAS *curfas = fas; 448656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 449656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 450656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 451656ede7eSPeter Brune bottom = y - th; 452656ede7eSPeter Brune while (curfas) { 453b4375e8dSPeter Brune if (!curfas->smoothu) { 454656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 455656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 456656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 457b4375e8dSPeter Brune } else { 458b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 459b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 460b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 461b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 462b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 463b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 464b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 465b4375e8dSPeter Brune } 466656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 467656ede7eSPeter Brune bottom -= 5*th; 4681aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 4690298fd71SBarry Smith else curfas = NULL; 470656ede7eSPeter Brune } 471421d9b32SPeter Brune } 472ab8d36c9SPeter Brune } 473421d9b32SPeter Brune PetscFunctionReturn(0); 474421d9b32SPeter Brune } 475421d9b32SPeter Brune 476421d9b32SPeter Brune #undef __FUNCT__ 47791f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private" 47839bd7f45SPeter Brune /* 47939bd7f45SPeter Brune Defines the action of the downsmoother 48039bd7f45SPeter Brune */ 48191f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 482b9c2fdf1SPeter Brune { 48339bd7f45SPeter Brune PetscErrorCode ierr = 0; 484742fe5e2SPeter Brune SNESConvergedReason reason; 485ab8d36c9SPeter Brune Vec FPC; 486ab8d36c9SPeter Brune SNES smoothd; 4870dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 4886e111a19SKarl Rupp 489421d9b32SPeter Brune PetscFunctionBegin; 490ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 491e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 4920dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 493ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 4940dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 495742fe5e2SPeter Brune /* check convergence reason for the smoother */ 496ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 497e70c42e5SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) { 498742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 499742fe5e2SPeter Brune PetscFunctionReturn(0); 500742fe5e2SPeter Brune } 5010298fd71SBarry Smith ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 5024b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 503b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr); 50439bd7f45SPeter Brune PetscFunctionReturn(0); 50539bd7f45SPeter Brune } 50639bd7f45SPeter Brune 50739bd7f45SPeter Brune 50839bd7f45SPeter Brune #undef __FUNCT__ 50991f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private" 51039bd7f45SPeter Brune /* 51107144faaSPeter Brune Defines the action of the upsmoother 51239bd7f45SPeter Brune */ 5130adebc6cSBarry Smith PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 5140adebc6cSBarry Smith { 51539bd7f45SPeter Brune PetscErrorCode ierr = 0; 51639bd7f45SPeter Brune SNESConvergedReason reason; 517ab8d36c9SPeter Brune Vec FPC; 518ab8d36c9SPeter Brune SNES smoothu; 5190dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 520ab8d36c9SPeter Brune 5216e111a19SKarl Rupp PetscFunctionBegin; 522ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 5230dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 524ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 5250dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 52639bd7f45SPeter Brune /* check convergence reason for the smoother */ 527ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 52839bd7f45SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 52939bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 53039bd7f45SPeter Brune PetscFunctionReturn(0); 53139bd7f45SPeter Brune } 5320298fd71SBarry Smith ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 5334b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 534b9c2fdf1SPeter Brune ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr); 53539bd7f45SPeter Brune PetscFunctionReturn(0); 53639bd7f45SPeter Brune } 53739bd7f45SPeter Brune 53839bd7f45SPeter Brune #undef __FUNCT__ 539938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec" 540938e4a01SJed Brown /*@ 541938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 542938e4a01SJed Brown 543938e4a01SJed Brown Collective 544938e4a01SJed Brown 545938e4a01SJed Brown Input Arguments: 546938e4a01SJed Brown . snes - SNESFAS 547938e4a01SJed Brown 548938e4a01SJed Brown Output Arguments: 549938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 550938e4a01SJed Brown 551938e4a01SJed Brown Level: developer 552938e4a01SJed Brown 553938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 554938e4a01SJed Brown @*/ 555938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 556938e4a01SJed Brown { 557938e4a01SJed Brown PetscErrorCode ierr; 558938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 559938e4a01SJed Brown 560938e4a01SJed Brown PetscFunctionBegin; 5611aa26658SKarl Rupp if (fas->rscale) { 5621aa26658SKarl Rupp ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 563f5af7f23SKarl Rupp } else if (fas->inject) { 5640298fd71SBarry Smith ierr = MatGetVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 565ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 566938e4a01SJed Brown PetscFunctionReturn(0); 567938e4a01SJed Brown } 568938e4a01SJed Brown 569e9923e8dSJed Brown #undef __FUNCT__ 570e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict" 571e9923e8dSJed Brown /*@ 572e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 573e9923e8dSJed Brown 574e9923e8dSJed Brown Collective 575e9923e8dSJed Brown 576e9923e8dSJed Brown Input Arguments: 577e9923e8dSJed Brown + fine - SNES from which to restrict 578e9923e8dSJed Brown - Xfine - vector to restrict 579e9923e8dSJed Brown 580e9923e8dSJed Brown Output Arguments: 581e9923e8dSJed Brown . Xcoarse - result of restriction 582e9923e8dSJed Brown 583e9923e8dSJed Brown Level: developer 584e9923e8dSJed Brown 585e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 586e9923e8dSJed Brown @*/ 587e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 588e9923e8dSJed Brown { 589e9923e8dSJed Brown PetscErrorCode ierr; 590e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 591e9923e8dSJed Brown 592e9923e8dSJed Brown PetscFunctionBegin; 593e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 594e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 595e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 596e9923e8dSJed Brown if (fas->inject) { 597e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 598e9923e8dSJed Brown } else { 599e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 600e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 601e9923e8dSJed Brown } 602e9923e8dSJed Brown PetscFunctionReturn(0); 603e9923e8dSJed Brown } 604e9923e8dSJed Brown 605e9923e8dSJed Brown #undef __FUNCT__ 6068c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection" 60739bd7f45SPeter Brune /* 60839bd7f45SPeter Brune 60939bd7f45SPeter Brune Performs the FAS coarse correction as: 61039bd7f45SPeter Brune 61139bd7f45SPeter Brune fine problem: F(x) = 0 61239bd7f45SPeter Brune coarse problem: F^c(x) = b^c 61339bd7f45SPeter Brune 61439bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 61539bd7f45SPeter Brune 61639bd7f45SPeter Brune */ 6170adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 6180adebc6cSBarry Smith { 61939bd7f45SPeter Brune PetscErrorCode ierr; 62039bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 62139bd7f45SPeter Brune SNESConvergedReason reason; 622ab8d36c9SPeter Brune SNES next; 623ab8d36c9SPeter Brune Mat restrct, interpolate; 6240dd27c6cSPeter Brune SNES_FAS *fasc; 6255fd66863SKarl Rupp 62639bd7f45SPeter Brune PetscFunctionBegin; 627ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 628ab8d36c9SPeter Brune if (next) { 6290dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 6300dd27c6cSPeter Brune 631ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 632ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 633ab8d36c9SPeter Brune 634ab8d36c9SPeter Brune X_c = next->vec_sol; 635ab8d36c9SPeter Brune Xo_c = next->work[0]; 636ab8d36c9SPeter Brune F_c = next->vec_func; 637ab8d36c9SPeter Brune B_c = next->vec_rhs; 638efe1f98aSPeter Brune 6390dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 640938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 641293a7e31SPeter Brune /* restrict the defect */ 642ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 6430dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 6440dd27c6cSPeter Brune 6450dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 646ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 6470dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 6480dd27c6cSPeter Brune 6490dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 650e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 651b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 652e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 653ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 654ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 655c90fad12SPeter Brune 656c90fad12SPeter Brune /* recurse to the next level */ 657e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 658ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 659ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 660742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 661742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 662742fe5e2SPeter Brune PetscFunctionReturn(0); 663742fe5e2SPeter Brune } 664fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 665fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 6660dd27c6cSPeter Brune 6670dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 668ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 6690dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 670293a7e31SPeter Brune } 67139bd7f45SPeter Brune PetscFunctionReturn(0); 67239bd7f45SPeter Brune } 67339bd7f45SPeter Brune 67439bd7f45SPeter Brune #undef __FUNCT__ 6752cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive" 67639bd7f45SPeter Brune /* 67739bd7f45SPeter Brune 67839bd7f45SPeter Brune The additive cycle looks like: 67939bd7f45SPeter Brune 68007144faaSPeter Brune xhat = x 68107144faaSPeter Brune xhat = dS(x, b) 68207144faaSPeter Brune x = coarsecorrection(xhat, b_d) 68307144faaSPeter Brune x = x + nu*(xhat - x); 68439bd7f45SPeter Brune (optional) x = uS(x, b) 68539bd7f45SPeter Brune 68639bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 68739bd7f45SPeter Brune 68839bd7f45SPeter Brune */ 6890adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 6900adebc6cSBarry Smith { 69107144faaSPeter Brune Vec F, B, Xhat; 69222c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 69339bd7f45SPeter Brune PetscErrorCode ierr; 69407144faaSPeter Brune SNESConvergedReason reason; 69522c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 69622c1e704SPeter Brune PetscBool lssuccess; 697ab8d36c9SPeter Brune SNES next; 698ab8d36c9SPeter Brune Mat restrct, interpolate; 6990dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 7000dd27c6cSPeter Brune 70139bd7f45SPeter Brune PetscFunctionBegin; 702ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 70339bd7f45SPeter Brune F = snes->vec_func; 70439bd7f45SPeter Brune B = snes->vec_rhs; 705e7f468e7SPeter Brune Xhat = snes->work[1]; 70607144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 70707144faaSPeter Brune /* recurse first */ 708ab8d36c9SPeter Brune if (next) { 7090dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 710ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 711ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 7120dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 71307144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 7140dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 715c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 716ab8d36c9SPeter Brune X_c = next->vec_sol; 717ab8d36c9SPeter Brune Xo_c = next->work[0]; 718ab8d36c9SPeter Brune F_c = next->vec_func; 719ab8d36c9SPeter Brune B_c = next->vec_rhs; 72039bd7f45SPeter Brune 721938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 72207144faaSPeter Brune /* restrict the defect */ 723ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 72407144faaSPeter Brune 72507144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 7260dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 727ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 7280dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 729e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 730b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 731e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 73207144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 73307144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 73407144faaSPeter Brune 73507144faaSPeter Brune /* recurse */ 736e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 737ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 73807144faaSPeter Brune 73907144faaSPeter Brune /* smooth on this level */ 74091f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 74107144faaSPeter Brune 742ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 74307144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 74407144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 74507144faaSPeter Brune PetscFunctionReturn(0); 74607144faaSPeter Brune } 74707144faaSPeter Brune 74807144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 749c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 750ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 75107144faaSPeter Brune 752ddebd997SPeter Brune /* additive correction of the coarse direction*/ 753f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 754f1c6b773SPeter Brune ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr); 7559e764e56SPeter Brune if (!lssuccess) { 7569e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 7579e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 7589e764e56SPeter Brune PetscFunctionReturn(0); 7599e764e56SPeter Brune } 7609e764e56SPeter Brune } 761b9c2fdf1SPeter Brune ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 76207144faaSPeter Brune } else { 76391f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 76407144faaSPeter Brune } 76539bd7f45SPeter Brune PetscFunctionReturn(0); 76639bd7f45SPeter Brune } 76739bd7f45SPeter Brune 76839bd7f45SPeter Brune #undef __FUNCT__ 7692cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative" 77039bd7f45SPeter Brune /* 77139bd7f45SPeter Brune 77239bd7f45SPeter Brune Defines the FAS cycle as: 77339bd7f45SPeter Brune 77439bd7f45SPeter Brune fine problem: F(x) = 0 77539bd7f45SPeter Brune coarse problem: F^c(x) = b^c 77639bd7f45SPeter Brune 77739bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x)) 77839bd7f45SPeter Brune 77939bd7f45SPeter Brune correction: 78039bd7f45SPeter Brune 78139bd7f45SPeter Brune x = x + I(x^c - Rx) 78239bd7f45SPeter Brune 78339bd7f45SPeter Brune */ 7840adebc6cSBarry Smith PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7850adebc6cSBarry Smith { 78639bd7f45SPeter Brune 78739bd7f45SPeter Brune PetscErrorCode ierr; 78839bd7f45SPeter Brune Vec F,B; 78939bd7f45SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 79039bd7f45SPeter Brune 79139bd7f45SPeter Brune PetscFunctionBegin; 79239bd7f45SPeter Brune F = snes->vec_func; 79339bd7f45SPeter Brune B = snes->vec_rhs; 79439bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 79591f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 796c90fad12SPeter Brune if (fas->level != 0) { 7978c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 79891f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 799fe6f9142SPeter Brune } 800fa9694d7SPeter Brune PetscFunctionReturn(0); 801421d9b32SPeter Brune } 802421d9b32SPeter Brune 803421d9b32SPeter Brune #undef __FUNCT__ 804421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS" 805421d9b32SPeter Brune 806421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes) 807421d9b32SPeter Brune { 808fa9694d7SPeter Brune PetscErrorCode ierr; 809fe6f9142SPeter Brune PetscInt i, maxits; 810ddb5aff1SPeter Brune Vec X, F; 811fe6f9142SPeter Brune PetscReal fnorm; 812b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 813b17ce1afSJed Brown DM dm; 814e70c42e5SPeter Brune PetscBool isFine; 815b17ce1afSJed Brown 816421d9b32SPeter Brune PetscFunctionBegin; 817fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 818fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 819fa9694d7SPeter Brune X = snes->vec_sol; 820f5a6d4f9SBarry Smith F = snes->vec_func; 821293a7e31SPeter Brune 82218a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 823293a7e31SPeter Brune /*norm setup */ 824ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 825fe6f9142SPeter Brune snes->iter = 0; 826fe6f9142SPeter Brune snes->norm = 0.; 827ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 828e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 8290dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 830fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 8310dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 832fe6f9142SPeter Brune if (snes->domainerror) { 833fe6f9142SPeter Brune snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN; 834fe6f9142SPeter Brune PetscFunctionReturn(0); 835fe6f9142SPeter Brune } 8361aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 837e4ed7901SPeter Brune 838fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 839189a9710SBarry Smith if (PetscIsInfOrNanReal(fnorm)) { 840189a9710SBarry Smith snes->reason = SNES_DIVERGED_FNORM_NAN; 841189a9710SBarry Smith PetscFunctionReturn(0); 842189a9710SBarry Smith } 843e4ed7901SPeter Brune 8446c4a5ce2SBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 845fe6f9142SPeter Brune snes->norm = fnorm; 846ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 847a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 848fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 849fe6f9142SPeter Brune 850fe6f9142SPeter Brune /* test convergence */ 851fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 852fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 853e4ed7901SPeter Brune 854b17ce1afSJed Brown 855b9c2fdf1SPeter Brune if (isFine) { 856b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 857b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 858b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 859b17ce1afSJed Brown DM dmcoarse; 860b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 861b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 862b17ce1afSJed Brown dm = dmcoarse; 863b17ce1afSJed Brown } 864b9c2fdf1SPeter Brune } 865b17ce1afSJed Brown 866fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 867fe6f9142SPeter Brune /* Call general purpose update function */ 868646217ecSPeter Brune 869fe6f9142SPeter Brune if (snes->ops->update) { 870fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 871fe6f9142SPeter Brune } 87207144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 87391f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 87407144faaSPeter Brune } else { 87591f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 87607144faaSPeter Brune } 877742fe5e2SPeter Brune 878742fe5e2SPeter Brune /* check for FAS cycle divergence */ 8791aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 880b9c2fdf1SPeter Brune 881c90fad12SPeter Brune /* Monitor convergence */ 882ce8c27fbSBarry Smith ierr = PetscObjectAMSTakeAccess((PetscObject)snes);CHKERRQ(ierr); 883c90fad12SPeter Brune snes->iter = i+1; 884ce8c27fbSBarry Smith ierr = PetscObjectAMSGrantAccess((PetscObject)snes);CHKERRQ(ierr); 885a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 886c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 887c90fad12SPeter Brune /* Test for convergence */ 88866585501SPeter Brune if (isFine) { 889b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 890c90fad12SPeter Brune if (snes->reason) break; 891fe6f9142SPeter Brune } 89266585501SPeter Brune } 893fe6f9142SPeter Brune if (i == maxits) { 894fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 895fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 896fe6f9142SPeter Brune } 897421d9b32SPeter Brune PetscFunctionReturn(0); 898421d9b32SPeter Brune } 899