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 66273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES,Vec,Vec,void*); 7efe1f98aSPeter Brune 840244768SBarry Smith static PetscErrorCode SNESReset_FAS(SNES snes) 9421d9b32SPeter Brune { 1077df8cc4SPeter Brune PetscErrorCode ierr = 0; 11421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 12421d9b32SPeter Brune 13421d9b32SPeter Brune PetscFunctionBegin; 14ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr); 15ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 163dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 17bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 18bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 19bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 207cd33a7bSPeter Brune if (fas->galerkin) { 217cd33a7bSPeter Brune ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr); 227cd33a7bSPeter Brune ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr); 237cd33a7bSPeter Brune } 24d477d801SPeter Brune if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);} 25421d9b32SPeter Brune PetscFunctionReturn(0); 26421d9b32SPeter Brune } 27421d9b32SPeter Brune 2840244768SBarry Smith static PetscErrorCode SNESDestroy_FAS(SNES snes) 29421d9b32SPeter Brune { 30421d9b32SPeter Brune SNES_FAS * fas = (SNES_FAS*)snes->data; 31742fe5e2SPeter Brune PetscErrorCode ierr = 0; 32421d9b32SPeter Brune 33421d9b32SPeter Brune PetscFunctionBegin; 34421d9b32SPeter Brune /* recursively resets and then destroys */ 3579d9a41aSPeter Brune ierr = SNESReset(snes);CHKERRQ(ierr); 361aa26658SKarl Rupp if (fas->next) { 371aa26658SKarl Rupp ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 381aa26658SKarl Rupp } 39421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 40421d9b32SPeter Brune PetscFunctionReturn(0); 41421d9b32SPeter Brune } 42421d9b32SPeter Brune 4340244768SBarry Smith static PetscErrorCode SNESSetUp_FAS(SNES snes) 44421d9b32SPeter Brune { 4548bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 46421d9b32SPeter Brune PetscErrorCode ierr; 47d1adcc6fSPeter Brune PetscInt dm_levels; 483dccd265SPeter Brune Vec vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */ 49ab8d36c9SPeter Brune SNES next; 500a7266b2SPatrick Farrell PetscBool isFine, hasCreateRestriction; 51f89ba88eSPeter Brune SNESLineSearch linesearch; 52f89ba88eSPeter Brune SNESLineSearch slinesearch; 53f89ba88eSPeter Brune void *lsprectx,*lspostctx; 546b2b7091SBarry Smith PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 556b2b7091SBarry Smith PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 56eff52c0eSPeter Brune 576b2b7091SBarry Smith PetscFunctionBegin; 58ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 59ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 60d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 61d1adcc6fSPeter Brune dm_levels++; 62cc05f883SPeter Brune if (dm_levels > fas->levels) { 632e8ce248SJed Brown /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/ 643dccd265SPeter Brune vec_sol = snes->vec_sol; 653dccd265SPeter Brune vec_func = snes->vec_func; 663dccd265SPeter Brune vec_sol_update = snes->vec_sol_update; 673dccd265SPeter Brune vec_rhs = snes->vec_rhs; 680298fd71SBarry Smith snes->vec_sol = NULL; 690298fd71SBarry Smith snes->vec_func = NULL; 700298fd71SBarry Smith snes->vec_sol_update = NULL; 710298fd71SBarry Smith snes->vec_rhs = NULL; 723dccd265SPeter Brune 733dccd265SPeter Brune /* reset the number of levels */ 740298fd71SBarry Smith ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr); 75cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 763dccd265SPeter Brune 773dccd265SPeter Brune snes->vec_sol = vec_sol; 783dccd265SPeter Brune snes->vec_func = vec_func; 793dccd265SPeter Brune snes->vec_rhs = vec_rhs; 803dccd265SPeter Brune snes->vec_sol_update = vec_sol_update; 81d1adcc6fSPeter Brune } 82d1adcc6fSPeter Brune } 83ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 84ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 853dccd265SPeter Brune 86fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 87cc05f883SPeter Brune 88ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 89ab8d36c9SPeter Brune if (!fas->smoothd) { 90ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 91ab8d36c9SPeter Brune } 92ab8d36c9SPeter Brune 9379d9a41aSPeter Brune if (snes->dm) { 94ab8d36c9SPeter Brune /* set the smoother DMs properly */ 95ab8d36c9SPeter Brune if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr); 96ab8d36c9SPeter Brune ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 9779d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 98ab8d36c9SPeter Brune if (next) { 9979d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 100ab8d36c9SPeter Brune if (!next->dm) { 101ce94432eSBarry Smith ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr); 102ab8d36c9SPeter Brune ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 10379d9a41aSPeter Brune } 10479d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 10579d9a41aSPeter Brune if (!fas->interpolate) { 106ab8d36c9SPeter Brune ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 107bccf9bb3SJed Brown if (!fas->restrct) { 1080a7266b2SPatrick Farrell ierr = DMHasCreateRestriction(next->dm, &hasCreateRestriction);CHKERRQ(ierr); 1090a7266b2SPatrick Farrell /* DM can create restrictions, use that */ 1100a7266b2SPatrick Farrell if (hasCreateRestriction) { 1110a7266b2SPatrick Farrell ierr = DMCreateRestriction(next->dm, snes->dm, &fas->restrct);CHKERRQ(ierr); 1120a7266b2SPatrick Farrell } else { 113bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 11479d9a41aSPeter Brune fas->restrct = fas->interpolate; 11579d9a41aSPeter Brune } 116bccf9bb3SJed Brown } 1170a7266b2SPatrick Farrell } 11879d9a41aSPeter Brune /* set the injection from the DM */ 11979d9a41aSPeter Brune if (!fas->inject) { 1206dbf9973SLawrence Mitchell ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr); 12179d9a41aSPeter Brune } 12279d9a41aSPeter Brune } 12379d9a41aSPeter Brune } 12479d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 12579d9a41aSPeter Brune if (fas->galerkin) { 1261aa26658SKarl Rupp if (next) { 1270298fd71SBarry Smith ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr); 1281aa26658SKarl Rupp } 1291aa26658SKarl Rupp if (fas->smoothd && fas->level != fas->levels - 1) { 1300298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 1311aa26658SKarl Rupp } 1321aa26658SKarl Rupp if (fas->smoothu && fas->level != fas->levels - 1) { 1330298fd71SBarry Smith ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr); 1341aa26658SKarl Rupp } 13579d9a41aSPeter Brune } 13679d9a41aSPeter Brune 137534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 138534ebe21SPeter Brune if (fas->smoothd) { 139bc3f2f05SPeter Brune if (fas->level == 0 && fas->levels != 1) { 140365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); 141534ebe21SPeter Brune } else { 142365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 143534ebe21SPeter Brune } 1447fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr); 145534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr); 1467601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 1477601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr); 1486b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 1496b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 1506b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 1516b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 152f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 1530dd27c6cSPeter Brune 1540dd27c6cSPeter Brune fas->smoothd->vec_sol = snes->vec_sol; 1550dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 1560dd27c6cSPeter Brune fas->smoothd->vec_sol_update = snes->vec_sol_update; 1570dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 1580dd27c6cSPeter Brune fas->smoothd->vec_func = snes->vec_func; 1590dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 1600dd27c6cSPeter Brune 1610dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1620dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr); 1630dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 164534ebe21SPeter Brune } 165534ebe21SPeter Brune 166534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 167534ebe21SPeter Brune if (fas->smoothu) { 168534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 169365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); 170534ebe21SPeter Brune } else { 171365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 172534ebe21SPeter Brune } 1737fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr); 174534ebe21SPeter Brune ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr); 1757601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 1767601faf0SJed Brown ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr); 1776b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 1786b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 1796b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 1806b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 181f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 1820dd27c6cSPeter Brune 1830dd27c6cSPeter Brune fas->smoothu->vec_sol = snes->vec_sol; 1840dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 1850dd27c6cSPeter Brune fas->smoothu->vec_sol_update = snes->vec_sol_update; 1860dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 1870dd27c6cSPeter Brune fas->smoothu->vec_func = snes->vec_func; 1880dd27c6cSPeter Brune ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 1890dd27c6cSPeter Brune 1900dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1910dd27c6cSPeter Brune ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr); 1920dd27c6cSPeter Brune if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);} 1930dd27c6cSPeter Brune 194534ebe21SPeter Brune } 195d06165b7SPeter Brune 196ab8d36c9SPeter Brune if (next) { 19779d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 198ab8d36c9SPeter Brune if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 199ab8d36c9SPeter Brune if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 2007fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); 2017601faf0SJed Brown ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 2027601faf0SJed Brown ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr); 2036b2b7091SBarry Smith ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 2046b2b7091SBarry Smith ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 2056b2b7091SBarry Smith ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 2066b2b7091SBarry Smith ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 207f89ba88eSPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 208ab8d36c9SPeter Brune ierr = SNESSetUp(next);CHKERRQ(ierr); 20979d9a41aSPeter Brune } 2106273346dSPeter Brune /* setup FAS work vectors */ 2116273346dSPeter Brune if (fas->galerkin) { 2126273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 2136273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 2146273346dSPeter Brune } 215421d9b32SPeter Brune PetscFunctionReturn(0); 216421d9b32SPeter Brune } 217421d9b32SPeter Brune 21840244768SBarry Smith static PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes) 219421d9b32SPeter Brune { 220ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 221ee78dd50SPeter Brune PetscInt levels = 1; 22287f44e3fSPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE; 223421d9b32SPeter Brune PetscErrorCode ierr; 22407144faaSPeter Brune SNESFASType fastype; 225fde0ff24SPeter Brune const char *optionsprefix; 226f1c6b773SPeter Brune SNESLineSearch linesearch; 22766585501SPeter Brune PetscInt m, n_up, n_down; 228ab8d36c9SPeter Brune SNES next; 229ab8d36c9SPeter Brune PetscBool isFine; 230421d9b32SPeter Brune 231421d9b32SPeter Brune PetscFunctionBegin; 232ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 233e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNESFAS Options-----------------------------------");CHKERRQ(ierr); 234ee78dd50SPeter Brune 235ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 236ab8d36c9SPeter Brune if (isFine) { 237ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 238c732cbdbSBarry Smith if (!flg && snes->dm) { 239c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 240c732cbdbSBarry Smith levels++; 241d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 242c732cbdbSBarry Smith } 2430298fd71SBarry Smith ierr = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr); 24407144faaSPeter Brune fastype = fas->fastype; 24507144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 24607144faaSPeter Brune if (flg) { 24707144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 24807144faaSPeter Brune } 249ee78dd50SPeter Brune 250fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 251ab8d36c9SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 252ab8d36c9SPeter Brune if (flg) { 253ab8d36c9SPeter Brune ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 254fde0ff24SPeter Brune } 25587f44e3fSPeter Brune ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr); 25687f44e3fSPeter Brune if (flg) { 25787f44e3fSPeter Brune ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr); 25887f44e3fSPeter Brune } 259fde0ff24SPeter Brune 260ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 261ab8d36c9SPeter Brune if (flg) { 262ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 263ab8d36c9SPeter Brune } 264ee78dd50SPeter Brune 265928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 266928e959bSPeter 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); 267928e959bSPeter Brune if (flg) {SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);} 268928e959bSPeter Brune } 269928e959bSPeter Brune 27066585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 271162d76ddSPeter Brune 27266585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 273162d76ddSPeter Brune 274d142ab34SLawrence Mitchell { 275d142ab34SLawrence Mitchell PetscViewer viewer; 276d142ab34SLawrence Mitchell PetscViewerFormat format; 277d142ab34SLawrence Mitchell ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix, 278d142ab34SLawrence Mitchell "-snes_fas_monitor",&viewer,&format,&monflg);CHKERRQ(ierr); 279d142ab34SLawrence Mitchell if (monflg) { 280d142ab34SLawrence Mitchell PetscViewerAndFormat *vf; 281d142ab34SLawrence Mitchell ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 282d142ab34SLawrence Mitchell ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 283d142ab34SLawrence Mitchell ierr = SNESFASSetMonitor(snes,vf,PETSC_TRUE);CHKERRQ(ierr); 284d142ab34SLawrence Mitchell } 285d142ab34SLawrence Mitchell } 2860dd27c6cSPeter Brune flg = PETSC_FALSE; 2870dd27c6cSPeter Brune monflg = PETSC_TRUE; 2880dd27c6cSPeter Brune ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 2890dd27c6cSPeter Brune if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 290ab8d36c9SPeter Brune } 291ee78dd50SPeter Brune 292421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 2938cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 294162d76ddSPeter Brune if (upflg) { 29566585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 296162d76ddSPeter Brune } 297162d76ddSPeter Brune if (downflg) { 29866585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 299162d76ddSPeter Brune } 300eff52c0eSPeter Brune 3019e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 3029e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 3039e764e56SPeter Brune if (!snes->linesearch) { 3047601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 3051a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 3069e764e56SPeter Brune } 3079e764e56SPeter Brune } 3089e764e56SPeter Brune 309ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 310ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 311ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 312421d9b32SPeter Brune PetscFunctionReturn(0); 313421d9b32SPeter Brune } 314421d9b32SPeter Brune 3159804daf3SBarry Smith #include <petscdraw.h> 31640244768SBarry Smith static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 317421d9b32SPeter Brune { 318421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 319656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 320ab8d36c9SPeter Brune PetscInt i; 321421d9b32SPeter Brune PetscErrorCode ierr; 322ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 323421d9b32SPeter Brune 324421d9b32SPeter Brune PetscFunctionBegin; 325ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 326ab8d36c9SPeter Brune if (isFine) { 327251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 328656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 329421d9b32SPeter Brune if (iascii) { 330*efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 331ab8d36c9SPeter Brune if (fas->galerkin) { 332ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 333421d9b32SPeter Brune } else { 334ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 335421d9b32SPeter Brune } 336ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 337ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 338ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 339ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 340ab8d36c9SPeter Brune if (!i) { 341ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 342421d9b32SPeter Brune } else { 343ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 344421d9b32SPeter Brune } 345ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 346166b3ea4SJed Brown if (smoothd) { 347ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 348166b3ea4SJed Brown } else { 349166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 350166b3ea4SJed Brown } 351ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 352ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 353ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 354ab8d36c9SPeter Brune } else if (i) { 355ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 356ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 357166b3ea4SJed Brown if (smoothu) { 358ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 359166b3ea4SJed Brown } else { 360166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 361166b3ea4SJed Brown } 362ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 363ab8d36c9SPeter Brune } 364ab8d36c9SPeter Brune } 365656ede7eSPeter Brune } else if (isdraw) { 366656ede7eSPeter Brune PetscDraw draw; 367b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 368656ede7eSPeter Brune SNES_FAS *curfas = fas; 369656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 370656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 371656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 372656ede7eSPeter Brune bottom = y - th; 373656ede7eSPeter Brune while (curfas) { 374b4375e8dSPeter Brune if (!curfas->smoothu) { 375656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 376656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 377656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 378b4375e8dSPeter Brune } else { 379b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 380b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 381b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 382b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 383b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 384b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 385b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 386b4375e8dSPeter Brune } 387656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 388656ede7eSPeter Brune bottom -= 5*th; 3891aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 3900298fd71SBarry Smith else curfas = NULL; 391656ede7eSPeter Brune } 392421d9b32SPeter Brune } 393ab8d36c9SPeter Brune } 394421d9b32SPeter Brune PetscFunctionReturn(0); 395421d9b32SPeter Brune } 396421d9b32SPeter Brune 39739bd7f45SPeter Brune /* 39839bd7f45SPeter Brune Defines the action of the downsmoother 39939bd7f45SPeter Brune */ 40040244768SBarry Smith static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 401b9c2fdf1SPeter Brune { 40239bd7f45SPeter Brune PetscErrorCode ierr = 0; 403742fe5e2SPeter Brune SNESConvergedReason reason; 404ab8d36c9SPeter Brune Vec FPC; 405ab8d36c9SPeter Brune SNES smoothd; 4066cbb2f26SLawrence Mitchell PetscBool flg; 4070dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 4086e111a19SKarl Rupp 409421d9b32SPeter Brune PetscFunctionBegin; 410ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 411e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 4120dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 413ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 4140dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 415742fe5e2SPeter Brune /* check convergence reason for the smoother */ 416ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 4171ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 418742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 419742fe5e2SPeter Brune PetscFunctionReturn(0); 420742fe5e2SPeter Brune } 4216cbb2f26SLawrence Mitchell 4220298fd71SBarry Smith ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 4236cbb2f26SLawrence Mitchell ierr = SNESGetAlwaysComputesFinalResidual(smoothd, &flg);CHKERRQ(ierr); 4246cbb2f26SLawrence Mitchell if (!flg) { 4256cbb2f26SLawrence Mitchell ierr = SNESComputeFunction(smoothd, X, FPC);CHKERRQ(ierr); 4266cbb2f26SLawrence Mitchell } 4274b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 42871dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 42939bd7f45SPeter Brune PetscFunctionReturn(0); 43039bd7f45SPeter Brune } 43139bd7f45SPeter Brune 43239bd7f45SPeter Brune 43339bd7f45SPeter Brune /* 43407144faaSPeter Brune Defines the action of the upsmoother 43539bd7f45SPeter Brune */ 43640244768SBarry Smith static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 4370adebc6cSBarry Smith { 43839bd7f45SPeter Brune PetscErrorCode ierr = 0; 43939bd7f45SPeter Brune SNESConvergedReason reason; 440ab8d36c9SPeter Brune Vec FPC; 441ab8d36c9SPeter Brune SNES smoothu; 4426cbb2f26SLawrence Mitchell PetscBool flg; 4430dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 444ab8d36c9SPeter Brune 4456e111a19SKarl Rupp PetscFunctionBegin; 446ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 4470dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 448ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 4490dd27c6cSPeter Brune if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);} 45039bd7f45SPeter Brune /* check convergence reason for the smoother */ 451ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 4521ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 45339bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 45439bd7f45SPeter Brune PetscFunctionReturn(0); 45539bd7f45SPeter Brune } 4560298fd71SBarry Smith ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 4576cbb2f26SLawrence Mitchell ierr = SNESGetAlwaysComputesFinalResidual(smoothu, &flg);CHKERRQ(ierr); 4586cbb2f26SLawrence Mitchell if (!flg) { 4596cbb2f26SLawrence Mitchell ierr = SNESComputeFunction(smoothu, X, FPC);CHKERRQ(ierr); 4606cbb2f26SLawrence Mitchell } 4614b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 46271dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 46339bd7f45SPeter Brune PetscFunctionReturn(0); 46439bd7f45SPeter Brune } 46539bd7f45SPeter Brune 466938e4a01SJed Brown /*@ 467938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 468938e4a01SJed Brown 469938e4a01SJed Brown Collective 470938e4a01SJed Brown 471938e4a01SJed Brown Input Arguments: 472938e4a01SJed Brown . snes - SNESFAS 473938e4a01SJed Brown 474938e4a01SJed Brown Output Arguments: 475938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 476938e4a01SJed Brown 477938e4a01SJed Brown Level: developer 478938e4a01SJed Brown 479938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 480938e4a01SJed Brown @*/ 481938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 482938e4a01SJed Brown { 483938e4a01SJed Brown PetscErrorCode ierr; 484938e4a01SJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data; 485938e4a01SJed Brown 486938e4a01SJed Brown PetscFunctionBegin; 4871aa26658SKarl Rupp if (fas->rscale) { 4881aa26658SKarl Rupp ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 489f5af7f23SKarl Rupp } else if (fas->inject) { 4902a7a6963SBarry Smith ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 491ce94432eSBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr); 492938e4a01SJed Brown PetscFunctionReturn(0); 493938e4a01SJed Brown } 494938e4a01SJed Brown 495e9923e8dSJed Brown /*@ 496e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 497e9923e8dSJed Brown 498e9923e8dSJed Brown Collective 499e9923e8dSJed Brown 500e9923e8dSJed Brown Input Arguments: 501e9923e8dSJed Brown + fine - SNES from which to restrict 502e9923e8dSJed Brown - Xfine - vector to restrict 503e9923e8dSJed Brown 504e9923e8dSJed Brown Output Arguments: 505e9923e8dSJed Brown . Xcoarse - result of restriction 506e9923e8dSJed Brown 507e9923e8dSJed Brown Level: developer 508e9923e8dSJed Brown 509e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 510e9923e8dSJed Brown @*/ 511e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 512e9923e8dSJed Brown { 513e9923e8dSJed Brown PetscErrorCode ierr; 514e9923e8dSJed Brown SNES_FAS *fas = (SNES_FAS*)fine->data; 515e9923e8dSJed Brown 516e9923e8dSJed Brown PetscFunctionBegin; 517e9923e8dSJed Brown PetscValidHeaderSpecific(fine,SNES_CLASSID,1); 518e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 519e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 520e9923e8dSJed Brown if (fas->inject) { 521e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 522e9923e8dSJed Brown } else { 523e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 524e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 525e9923e8dSJed Brown } 526e9923e8dSJed Brown PetscFunctionReturn(0); 527e9923e8dSJed Brown } 528e9923e8dSJed Brown 52939bd7f45SPeter Brune /* 53039bd7f45SPeter Brune 53139bd7f45SPeter Brune Performs the FAS coarse correction as: 53239bd7f45SPeter Brune 533b5527d98SMatthew G. Knepley fine problem: F(x) = b 534b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c 53539bd7f45SPeter Brune 536b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 53739bd7f45SPeter Brune 53839bd7f45SPeter Brune */ 5390adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 5400adebc6cSBarry Smith { 54139bd7f45SPeter Brune PetscErrorCode ierr; 54239bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 54339bd7f45SPeter Brune SNESConvergedReason reason; 544ab8d36c9SPeter Brune SNES next; 545ab8d36c9SPeter Brune Mat restrct, interpolate; 5460dd27c6cSPeter Brune SNES_FAS *fasc; 5475fd66863SKarl Rupp 54839bd7f45SPeter Brune PetscFunctionBegin; 549ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 550ab8d36c9SPeter Brune if (next) { 5510dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 5520dd27c6cSPeter Brune 553ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 554ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 555ab8d36c9SPeter Brune 556ab8d36c9SPeter Brune X_c = next->vec_sol; 557ab8d36c9SPeter Brune Xo_c = next->work[0]; 558ab8d36c9SPeter Brune F_c = next->vec_func; 559ab8d36c9SPeter Brune B_c = next->vec_rhs; 560efe1f98aSPeter Brune 5610dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 562938e4a01SJed Brown ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr); 5635609cbf2SMatthew G. Knepley /* restrict the defect: R(F(x) - b) */ 564ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 5650dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 5660dd27c6cSPeter Brune 5670dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 5685609cbf2SMatthew G. Knepley /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */ 569ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 5700dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 5710dd27c6cSPeter Brune 5720dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 573e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 574b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 575e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 576ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 577ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 578c90fad12SPeter Brune 579c90fad12SPeter Brune /* recurse to the next level */ 580e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 581ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 582ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 583742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 584742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 585742fe5e2SPeter Brune PetscFunctionReturn(0); 586742fe5e2SPeter Brune } 587fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 588fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 5890dd27c6cSPeter Brune 5900dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 591ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 5920dd27c6cSPeter Brune if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);} 593293a7e31SPeter Brune } 59439bd7f45SPeter Brune PetscFunctionReturn(0); 59539bd7f45SPeter Brune } 59639bd7f45SPeter Brune 59739bd7f45SPeter Brune /* 59839bd7f45SPeter Brune 59939bd7f45SPeter Brune The additive cycle looks like: 60039bd7f45SPeter Brune 60107144faaSPeter Brune xhat = x 60207144faaSPeter Brune xhat = dS(x, b) 60307144faaSPeter Brune x = coarsecorrection(xhat, b_d) 60407144faaSPeter Brune x = x + nu*(xhat - x); 60539bd7f45SPeter Brune (optional) x = uS(x, b) 60639bd7f45SPeter Brune 60739bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 60839bd7f45SPeter Brune 60939bd7f45SPeter Brune */ 61040244768SBarry Smith static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 6110adebc6cSBarry Smith { 61207144faaSPeter Brune Vec F, B, Xhat; 61322c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 61439bd7f45SPeter Brune PetscErrorCode ierr; 61507144faaSPeter Brune SNESConvergedReason reason; 61622c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 617422a814eSBarry Smith SNESLineSearchReason lsresult; 618ab8d36c9SPeter Brune SNES next; 619ab8d36c9SPeter Brune Mat restrct, interpolate; 6200dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 6210dd27c6cSPeter Brune 62239bd7f45SPeter Brune PetscFunctionBegin; 623ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 62439bd7f45SPeter Brune F = snes->vec_func; 62539bd7f45SPeter Brune B = snes->vec_rhs; 626e7f468e7SPeter Brune Xhat = snes->work[1]; 62707144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 62807144faaSPeter Brune /* recurse first */ 629ab8d36c9SPeter Brune if (next) { 6300dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 631ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 632ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 6330dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 63407144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 6350dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 636c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 637ab8d36c9SPeter Brune X_c = next->vec_sol; 638ab8d36c9SPeter Brune Xo_c = next->work[0]; 639ab8d36c9SPeter Brune F_c = next->vec_func; 640ab8d36c9SPeter Brune B_c = next->vec_rhs; 64139bd7f45SPeter Brune 642938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 64307144faaSPeter Brune /* restrict the defect */ 644ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 64507144faaSPeter Brune 64607144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 6470dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 648ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 6490dd27c6cSPeter Brune if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);} 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); 65307144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 65407144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 65507144faaSPeter Brune 65607144faaSPeter Brune /* recurse */ 657e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 658ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 65907144faaSPeter Brune 66007144faaSPeter Brune /* smooth on this level */ 66191f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 66207144faaSPeter Brune 663ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 66407144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 66507144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 66607144faaSPeter Brune PetscFunctionReturn(0); 66707144faaSPeter Brune } 66807144faaSPeter Brune 66907144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 670c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 671ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 67207144faaSPeter Brune 673ddebd997SPeter Brune /* additive correction of the coarse direction*/ 674f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 675422a814eSBarry Smith ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr); 676422a814eSBarry Smith ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 677422a814eSBarry Smith if (lsresult) { 6789e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 6799e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 6809e764e56SPeter Brune PetscFunctionReturn(0); 6819e764e56SPeter Brune } 6829e764e56SPeter Brune } 68307144faaSPeter Brune } else { 68491f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 68507144faaSPeter Brune } 68639bd7f45SPeter Brune PetscFunctionReturn(0); 68739bd7f45SPeter Brune } 68839bd7f45SPeter Brune 68939bd7f45SPeter Brune /* 69039bd7f45SPeter Brune 69139bd7f45SPeter Brune Defines the FAS cycle as: 69239bd7f45SPeter Brune 693b5527d98SMatthew G. Knepley fine problem: F(x) = b 69439bd7f45SPeter Brune coarse problem: F^c(x) = b^c 69539bd7f45SPeter Brune 696b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 69739bd7f45SPeter Brune 69839bd7f45SPeter Brune correction: 69939bd7f45SPeter Brune 70039bd7f45SPeter Brune x = x + I(x^c - Rx) 70139bd7f45SPeter Brune 70239bd7f45SPeter Brune */ 70340244768SBarry Smith static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7040adebc6cSBarry Smith { 70539bd7f45SPeter Brune 70639bd7f45SPeter Brune PetscErrorCode ierr; 70739bd7f45SPeter Brune Vec F,B; 70834d65b3cSPeter Brune SNES next; 70939bd7f45SPeter Brune 71039bd7f45SPeter Brune PetscFunctionBegin; 71139bd7f45SPeter Brune F = snes->vec_func; 71239bd7f45SPeter Brune B = snes->vec_rhs; 71339bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 71434d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 71591f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 71634d65b3cSPeter Brune if (next) { 7178c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 71891f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 719fe6f9142SPeter Brune } 720fa9694d7SPeter Brune PetscFunctionReturn(0); 721421d9b32SPeter Brune } 722421d9b32SPeter Brune 72340244768SBarry Smith static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 7248c94862eSPeter Brune { 7258c94862eSPeter Brune SNES next; 7268c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 7278c94862eSPeter Brune PetscBool isFine; 7288c94862eSPeter Brune PetscErrorCode ierr; 7298c94862eSPeter Brune 7308c94862eSPeter Brune PetscFunctionBegin; 7318c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 7328c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 7338c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 7348c94862eSPeter Brune fas->full_stage = 0; 7358c94862eSPeter Brune if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);} 7368c94862eSPeter Brune PetscFunctionReturn(0); 7378c94862eSPeter Brune } 7388c94862eSPeter Brune 73940244768SBarry Smith static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 7408c94862eSPeter Brune { 7418c94862eSPeter Brune PetscErrorCode ierr; 7428c94862eSPeter Brune Vec F,B; 7438c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 7448c94862eSPeter Brune PetscBool isFine; 7458c94862eSPeter Brune SNES next; 7468c94862eSPeter Brune 7478c94862eSPeter Brune PetscFunctionBegin; 7488c94862eSPeter Brune F = snes->vec_func; 7498c94862eSPeter Brune B = snes->vec_rhs; 7508c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 7518c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 7528c94862eSPeter Brune 7538c94862eSPeter Brune if (isFine) { 7548c94862eSPeter Brune ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr); 7558c94862eSPeter Brune } 7568c94862eSPeter Brune 7578c94862eSPeter Brune if (fas->full_stage == 0) { 758928e959bSPeter Brune /* downsweep */ 7598c94862eSPeter Brune if (next) { 7608c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 761928e959bSPeter Brune if (fas->full_downsweep||isFine) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 7628c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 7638c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 7648c94862eSPeter Brune if (fas->level != 1) next->max_its -= 1; 7658c94862eSPeter Brune } else { 766a3a80b83SMatthew G. Knepley /* The smoother on the coarse level is the coarse solver */ 7678c94862eSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 7688c94862eSPeter Brune } 7698c94862eSPeter Brune fas->full_stage = 1; 7708c94862eSPeter Brune } else if (fas->full_stage == 1) { 7718c94862eSPeter Brune if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 7728c94862eSPeter Brune if (next) { 7738c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 7748c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 7758c94862eSPeter Brune } 7768c94862eSPeter Brune } 7778c94862eSPeter Brune /* final v-cycle */ 7788c94862eSPeter Brune if (isFine) { 7798c94862eSPeter Brune if (next) { 7808c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 7818c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 7828c94862eSPeter Brune } 7838c94862eSPeter Brune } 7848c94862eSPeter Brune PetscFunctionReturn(0); 7858c94862eSPeter Brune } 7868c94862eSPeter Brune 78740244768SBarry Smith static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 78834d65b3cSPeter Brune { 78934d65b3cSPeter Brune PetscErrorCode ierr; 79034d65b3cSPeter Brune Vec F,B; 79134d65b3cSPeter Brune SNES next; 79234d65b3cSPeter Brune 79334d65b3cSPeter Brune PetscFunctionBegin; 79434d65b3cSPeter Brune F = snes->vec_func; 79534d65b3cSPeter Brune B = snes->vec_rhs; 79634d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 79734d65b3cSPeter Brune if (next) { 79834d65b3cSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 79934d65b3cSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 80034d65b3cSPeter Brune } else { 80134d65b3cSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 80234d65b3cSPeter Brune } 80334d65b3cSPeter Brune PetscFunctionReturn(0); 80434d65b3cSPeter Brune } 80534d65b3cSPeter Brune 806fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE; 807fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n" 808fffbeea8SBarry Smith " title = {Composing Scalable Nonlinear Algebraic Solvers},\n" 809fffbeea8SBarry Smith " author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n" 810fffbeea8SBarry Smith " year = 2013,\n" 811fffbeea8SBarry Smith " type = Preprint,\n" 812fffbeea8SBarry Smith " number = {ANL/MCS-P2010-0112},\n" 813fffbeea8SBarry Smith " institution = {Argonne National Laboratory}\n}\n"; 814fffbeea8SBarry Smith 81540244768SBarry Smith static PetscErrorCode SNESSolve_FAS(SNES snes) 816421d9b32SPeter Brune { 817fa9694d7SPeter Brune PetscErrorCode ierr; 818fe6f9142SPeter Brune PetscInt i, maxits; 819ddb5aff1SPeter Brune Vec X, F; 820fe6f9142SPeter Brune PetscReal fnorm; 821b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 822b17ce1afSJed Brown DM dm; 823e70c42e5SPeter Brune PetscBool isFine; 824b17ce1afSJed Brown 825421d9b32SPeter Brune PetscFunctionBegin; 826c579b300SPatrick Farrell 8276c4ed002SBarry 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); 828c579b300SPatrick Farrell 829fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 830fe6f9142SPeter Brune maxits = snes->max_its; /* maximum number of iterations */ 831fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 832fa9694d7SPeter Brune X = snes->vec_sol; 833f5a6d4f9SBarry Smith F = snes->vec_func; 834293a7e31SPeter Brune 83518a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 836293a7e31SPeter Brune /*norm setup */ 837e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 838fe6f9142SPeter Brune snes->iter = 0; 839fe6f9142SPeter Brune snes->norm = 0.; 840e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 841e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 8420dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 843fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 8440dd27c6cSPeter Brune if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);} 8451aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 846e4ed7901SPeter Brune 847fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 848422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 849e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 850fe6f9142SPeter Brune snes->norm = fnorm; 851e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 852a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 853fe6f9142SPeter Brune ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr); 854fe6f9142SPeter Brune 855fe6f9142SPeter Brune /* test convergence */ 856fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 857fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 858e4ed7901SPeter Brune 859b17ce1afSJed Brown 860b9c2fdf1SPeter Brune if (isFine) { 861b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 862b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 863b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 864b17ce1afSJed Brown DM dmcoarse; 865b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 866b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 867b17ce1afSJed Brown dm = dmcoarse; 868b17ce1afSJed Brown } 869b9c2fdf1SPeter Brune } 870b17ce1afSJed Brown 871fe6f9142SPeter Brune for (i = 0; i < maxits; i++) { 872fe6f9142SPeter Brune /* Call general purpose update function */ 873646217ecSPeter Brune 874fe6f9142SPeter Brune if (snes->ops->update) { 875fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 876fe6f9142SPeter Brune } 87707144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 87891f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 8798c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_ADDITIVE) { 88091f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 8818c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_FULL) { 8828c94862eSPeter Brune ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr); 88334d65b3cSPeter Brune } else if (fas->fastype ==SNES_FAS_KASKADE) { 88434d65b3cSPeter Brune ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr); 8856c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type"); 886742fe5e2SPeter Brune 887742fe5e2SPeter Brune /* check for FAS cycle divergence */ 8881aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 889b9c2fdf1SPeter Brune 890c90fad12SPeter Brune /* Monitor convergence */ 891e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 892c90fad12SPeter Brune snes->iter = i+1; 893e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 894a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 895c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 896c90fad12SPeter Brune /* Test for convergence */ 89766585501SPeter Brune if (isFine) { 898b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 899c90fad12SPeter Brune if (snes->reason) break; 900fe6f9142SPeter Brune } 90166585501SPeter Brune } 902fe6f9142SPeter Brune if (i == maxits) { 903fe6f9142SPeter Brune ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr); 904fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 905fe6f9142SPeter Brune } 906421d9b32SPeter Brune PetscFunctionReturn(0); 907421d9b32SPeter Brune } 90840244768SBarry Smith 90940244768SBarry Smith /*MC 91040244768SBarry Smith 91140244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 91240244768SBarry Smith 91340244768SBarry Smith The nonlinear problem is solved by correction using coarse versions 91440244768SBarry Smith of the nonlinear problem. This problem is perturbed so that a projected 91540244768SBarry Smith solution of the fine problem elicits no correction from the coarse problem. 91640244768SBarry Smith 91740244768SBarry Smith Options Database: 91840244768SBarry Smith + -snes_fas_levels - The number of levels 91940244768SBarry Smith . -snes_fas_cycles<1> - The number of cycles -- 1 for V, 2 for W 92040244768SBarry Smith . -snes_fas_type<additive,multiplicative,full,kaskade> - Additive or multiplicative cycle 92140244768SBarry Smith . -snes_fas_galerkin<PETSC_FALSE> - Form coarse problems by projection back upon the fine problem 92240244768SBarry Smith . -snes_fas_smoothup<1> - The number of iterations of the post-smoother 92340244768SBarry Smith . -snes_fas_smoothdown<1> - The number of iterations of the pre-smoother 92440244768SBarry Smith . -snes_fas_monitor - Monitor progress of all of the levels 92540244768SBarry Smith . -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS 92640244768SBarry Smith . -fas_levels_snes_ - SNES options for all smoothers 92740244768SBarry Smith . -fas_levels_cycle_snes_ - SNES options for all cycles 92840244768SBarry Smith . -fas_levels_i_snes_ - SNES options for the smoothers on level i 92940244768SBarry Smith . -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i 93040244768SBarry Smith - -fas_coarse_snes_ - SNES options for the coarsest smoother 93140244768SBarry Smith 93240244768SBarry Smith Notes: 93340244768SBarry Smith The organization of the FAS solver is slightly different from the organization of PCMG 93440244768SBarry Smith As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance. 93540244768SBarry Smith The cycle SNES instance may be used for monitoring convergence on a particular level. 93640244768SBarry Smith 93740244768SBarry Smith Level: beginner 93840244768SBarry Smith 93940244768SBarry Smith References: 94040244768SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 94140244768SBarry Smith SIAM Review, 57(4), 2015 94240244768SBarry Smith 94340244768SBarry Smith .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 94440244768SBarry Smith M*/ 94540244768SBarry Smith 94640244768SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes) 94740244768SBarry Smith { 94840244768SBarry Smith SNES_FAS *fas; 94940244768SBarry Smith PetscErrorCode ierr; 95040244768SBarry Smith 95140244768SBarry Smith PetscFunctionBegin; 95240244768SBarry Smith snes->ops->destroy = SNESDestroy_FAS; 95340244768SBarry Smith snes->ops->setup = SNESSetUp_FAS; 95440244768SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_FAS; 95540244768SBarry Smith snes->ops->view = SNESView_FAS; 95640244768SBarry Smith snes->ops->solve = SNESSolve_FAS; 95740244768SBarry Smith snes->ops->reset = SNESReset_FAS; 95840244768SBarry Smith 95940244768SBarry Smith snes->usesksp = PETSC_FALSE; 960*efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 96140244768SBarry Smith 96240244768SBarry Smith if (!snes->tolerancesset) { 96340244768SBarry Smith snes->max_funcs = 30000; 96440244768SBarry Smith snes->max_its = 10000; 96540244768SBarry Smith } 96640244768SBarry Smith 9674fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 9684fc747eaSLawrence Mitchell 96940244768SBarry Smith ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr); 97040244768SBarry Smith 97140244768SBarry Smith snes->data = (void*) fas; 97240244768SBarry Smith fas->level = 0; 97340244768SBarry Smith fas->levels = 1; 97440244768SBarry Smith fas->n_cycles = 1; 97540244768SBarry Smith fas->max_up_it = 1; 97640244768SBarry Smith fas->max_down_it = 1; 97740244768SBarry Smith fas->smoothu = NULL; 97840244768SBarry Smith fas->smoothd = NULL; 97940244768SBarry Smith fas->next = NULL; 98040244768SBarry Smith fas->previous = NULL; 98140244768SBarry Smith fas->fine = snes; 98240244768SBarry Smith fas->interpolate = NULL; 98340244768SBarry Smith fas->restrct = NULL; 98440244768SBarry Smith fas->inject = NULL; 98540244768SBarry Smith fas->usedmfornumberoflevels = PETSC_FALSE; 98640244768SBarry Smith fas->fastype = SNES_FAS_MULTIPLICATIVE; 98740244768SBarry Smith fas->full_downsweep = PETSC_FALSE; 98840244768SBarry Smith 98940244768SBarry Smith fas->eventsmoothsetup = 0; 99040244768SBarry Smith fas->eventsmoothsolve = 0; 99140244768SBarry Smith fas->eventresidual = 0; 99240244768SBarry Smith fas->eventinterprestrict = 0; 99340244768SBarry Smith PetscFunctionReturn(0); 99440244768SBarry Smith } 995