1421d9b32SPeter Brune /* Defines the basic SNES object */ 2a055c8caSBarry Smith #include <../src/snes/impls/fas/fasimpls.h> /*I "petscsnes.h" I*/ 3421d9b32SPeter Brune 49e5d0892SLisandro Dalcin const char *const SNESFASTypes[] = {"MULTIPLICATIVE", "ADDITIVE", "FULL", "KASKADE", "SNESFASType", "SNES_FAS", NULL}; 507144faaSPeter Brune 69371c9d4SSatish Balay static PetscErrorCode SNESReset_FAS(SNES snes) { 7421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 8421d9b32SPeter Brune 9421d9b32SPeter Brune PetscFunctionBegin; 109566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&fas->smoothu)); 119566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&fas->smoothd)); 129566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->inject)); 139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->interpolate)); 149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->restrct)); 159566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->rscale)); 169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->Xg)); 179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->Fg)); 189566063dSJacob Faibussowitsch if (fas->next) PetscCall(SNESReset(fas->next)); 19421d9b32SPeter Brune PetscFunctionReturn(0); 20421d9b32SPeter Brune } 21421d9b32SPeter Brune 229371c9d4SSatish Balay static PetscErrorCode SNESDestroy_FAS(SNES snes) { 23421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 24421d9b32SPeter Brune 25421d9b32SPeter Brune PetscFunctionBegin; 26421d9b32SPeter Brune /* recursively resets and then destroys */ 279566063dSJacob Faibussowitsch PetscCall(SNESReset_FAS(snes)); 289566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&fas->next)); 299566063dSJacob Faibussowitsch PetscCall(PetscFree(fas)); 30421d9b32SPeter Brune PetscFunctionReturn(0); 31421d9b32SPeter Brune } 32421d9b32SPeter Brune 339371c9d4SSatish Balay static PetscErrorCode SNESFASSetUpLineSearch_Private(SNES snes, SNES smooth) { 34f833ba53SLisandro Dalcin SNESLineSearch linesearch; 35f833ba53SLisandro Dalcin SNESLineSearch slinesearch; 36f833ba53SLisandro Dalcin void *lsprectx, *lspostctx; 37f833ba53SLisandro Dalcin PetscErrorCode (*precheck)(SNESLineSearch, Vec, Vec, PetscBool *, void *); 38f833ba53SLisandro Dalcin PetscErrorCode (*postcheck)(SNESLineSearch, Vec, Vec, Vec, PetscBool *, PetscBool *, void *); 39f833ba53SLisandro Dalcin 40f833ba53SLisandro Dalcin PetscFunctionBegin; 41f833ba53SLisandro Dalcin if (!snes->linesearch) PetscFunctionReturn(0); 429566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 439566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(smooth, &slinesearch)); 449566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetPreCheck(linesearch, &precheck, &lsprectx)); 459566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetPostCheck(linesearch, &postcheck, &lspostctx)); 469566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(slinesearch, precheck, lsprectx)); 479566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPostCheck(slinesearch, postcheck, lspostctx)); 489566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch)); 49f833ba53SLisandro Dalcin PetscFunctionReturn(0); 50f833ba53SLisandro Dalcin } 51f833ba53SLisandro Dalcin 529371c9d4SSatish Balay static PetscErrorCode SNESFASCycleSetUpSmoother_Private(SNES snes, SNES smooth) { 53f833ba53SLisandro Dalcin SNES_FAS *fas = (SNES_FAS *)snes->data; 54f833ba53SLisandro Dalcin 55f833ba53SLisandro Dalcin PetscFunctionBegin; 569566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)smooth)); 579566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(smooth)); 589566063dSJacob Faibussowitsch PetscCall(SNESFASSetUpLineSearch_Private(snes, smooth)); 59f833ba53SLisandro Dalcin 609566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)snes->vec_sol)); 619566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)snes->vec_sol_update)); 629566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)snes->vec_func)); 63f833ba53SLisandro Dalcin smooth->vec_sol = snes->vec_sol; 64f833ba53SLisandro Dalcin smooth->vec_sol_update = snes->vec_sol_update; 65f833ba53SLisandro Dalcin smooth->vec_func = snes->vec_func; 66f833ba53SLisandro Dalcin 679566063dSJacob Faibussowitsch if (fas->eventsmoothsetup) PetscCall(PetscLogEventBegin(fas->eventsmoothsetup, smooth, 0, 0, 0)); 689566063dSJacob Faibussowitsch PetscCall(SNESSetUp(smooth)); 699566063dSJacob Faibussowitsch if (fas->eventsmoothsetup) PetscCall(PetscLogEventEnd(fas->eventsmoothsetup, smooth, 0, 0, 0)); 70f833ba53SLisandro Dalcin PetscFunctionReturn(0); 71f833ba53SLisandro Dalcin } 72f833ba53SLisandro Dalcin 739371c9d4SSatish Balay static PetscErrorCode SNESSetUp_FAS(SNES snes) { 7448bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 75d1adcc6fSPeter Brune PetscInt dm_levels; 76ab8d36c9SPeter Brune SNES next; 77f833ba53SLisandro Dalcin PetscBool isFine, hasCreateRestriction, hasCreateInjection; 78eff52c0eSPeter Brune 796b2b7091SBarry Smith PetscFunctionBegin; 809566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 81ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 829566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(snes->dm, &dm_levels)); 83d1adcc6fSPeter Brune dm_levels++; 84cc05f883SPeter Brune if (dm_levels > fas->levels) { 853dccd265SPeter Brune /* reset the number of levels */ 869566063dSJacob Faibussowitsch PetscCall(SNESFASSetLevels(snes, dm_levels, NULL)); 879566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 88d1adcc6fSPeter Brune } 89d1adcc6fSPeter Brune } 909566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 91ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 923dccd265SPeter Brune 939566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 2)); /* work vectors used for intergrid transfers */ 94cc05f883SPeter Brune 95ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 96*48a46eb9SPierre Jolivet if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd)); 97ab8d36c9SPeter Brune 9879d9a41aSPeter Brune if (snes->dm) { 99ab8d36c9SPeter Brune /* set the smoother DMs properly */ 1009566063dSJacob Faibussowitsch if (fas->smoothu) PetscCall(SNESSetDM(fas->smoothu, snes->dm)); 1019566063dSJacob Faibussowitsch PetscCall(SNESSetDM(fas->smoothd, snes->dm)); 10279d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 103ab8d36c9SPeter Brune if (next) { 10479d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 105ab8d36c9SPeter Brune if (!next->dm) { 1069566063dSJacob Faibussowitsch PetscCall(DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm)); 1079566063dSJacob Faibussowitsch PetscCall(SNESSetDM(next, next->dm)); 10879d9a41aSPeter Brune } 10979d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 11079d9a41aSPeter Brune if (!fas->interpolate) { 1119566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale)); 112bccf9bb3SJed Brown if (!fas->restrct) { 1139566063dSJacob Faibussowitsch PetscCall(DMHasCreateRestriction(next->dm, &hasCreateRestriction)); 1140a7266b2SPatrick Farrell /* DM can create restrictions, use that */ 1150a7266b2SPatrick Farrell if (hasCreateRestriction) { 1169566063dSJacob Faibussowitsch PetscCall(DMCreateRestriction(next->dm, snes->dm, &fas->restrct)); 1170a7266b2SPatrick Farrell } else { 1189566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fas->interpolate)); 11979d9a41aSPeter Brune fas->restrct = fas->interpolate; 12079d9a41aSPeter Brune } 121bccf9bb3SJed Brown } 1220a7266b2SPatrick Farrell } 12379d9a41aSPeter Brune /* set the injection from the DM */ 12479d9a41aSPeter Brune if (!fas->inject) { 1259566063dSJacob Faibussowitsch PetscCall(DMHasCreateInjection(next->dm, &hasCreateInjection)); 126*48a46eb9SPierre Jolivet if (hasCreateInjection) PetscCall(DMCreateInjection(next->dm, snes->dm, &fas->inject)); 12779d9a41aSPeter Brune } 12879d9a41aSPeter Brune } 12923e68893SLawrence Mitchell } 130f833ba53SLisandro Dalcin 13179d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 13279d9a41aSPeter Brune if (fas->galerkin) { 1331baa6e33SBarry Smith if (next) PetscCall(SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next)); 134*48a46eb9SPierre Jolivet if (fas->smoothd && fas->level != fas->levels - 1) PetscCall(SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes)); 135*48a46eb9SPierre Jolivet if (fas->smoothu && fas->level != fas->levels - 1) PetscCall(SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes)); 13679d9a41aSPeter Brune } 13779d9a41aSPeter Brune 138534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 139534ebe21SPeter Brune if (fas->smoothd) { 140bc3f2f05SPeter Brune if (fas->level == 0 && fas->levels != 1) { 1419566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE)); 142534ebe21SPeter Brune } else { 1439566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY)); 144534ebe21SPeter Brune } 1459566063dSJacob Faibussowitsch PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothd)); 146534ebe21SPeter Brune } 147534ebe21SPeter Brune 148534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 149534ebe21SPeter Brune if (fas->smoothu) { 150534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 1519566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE)); 152534ebe21SPeter Brune } else { 1539566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY)); 154534ebe21SPeter Brune } 1559566063dSJacob Faibussowitsch PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothu)); 156534ebe21SPeter Brune } 157d06165b7SPeter Brune 158ab8d36c9SPeter Brune if (next) { 15979d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 1609566063dSJacob Faibussowitsch if (!next->vec_sol) PetscCall(SNESFASCreateCoarseVec(snes, &next->vec_sol)); 1619566063dSJacob Faibussowitsch if (!next->vec_rhs) PetscCall(SNESFASCreateCoarseVec(snes, &next->vec_rhs)); 1629566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next)); 1639566063dSJacob Faibussowitsch PetscCall(SNESFASSetUpLineSearch_Private(snes, next)); 1649566063dSJacob Faibussowitsch PetscCall(SNESSetUp(next)); 16579d9a41aSPeter Brune } 166f833ba53SLisandro Dalcin 1676273346dSPeter Brune /* setup FAS work vectors */ 1686273346dSPeter Brune if (fas->galerkin) { 1699566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &fas->Xg)); 1709566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &fas->Fg)); 1716273346dSPeter Brune } 172421d9b32SPeter Brune PetscFunctionReturn(0); 173421d9b32SPeter Brune } 174421d9b32SPeter Brune 1759371c9d4SSatish Balay static PetscErrorCode SNESSetFromOptions_FAS(SNES snes, PetscOptionItems *PetscOptionsObject) { 176ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 177ee78dd50SPeter Brune PetscInt levels = 1; 17887f44e3fSPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE, continuationflg = PETSC_FALSE; 17907144faaSPeter Brune SNESFASType fastype; 180fde0ff24SPeter Brune const char *optionsprefix; 181f1c6b773SPeter Brune SNESLineSearch linesearch; 18266585501SPeter Brune PetscInt m, n_up, n_down; 183ab8d36c9SPeter Brune SNES next; 184ab8d36c9SPeter Brune PetscBool isFine; 185421d9b32SPeter Brune 186421d9b32SPeter Brune PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 188d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject, "SNESFAS Options-----------------------------------"); 189ee78dd50SPeter Brune 190ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 191ab8d36c9SPeter Brune if (isFine) { 1929566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg)); 193c732cbdbSBarry Smith if (!flg && snes->dm) { 1949566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(snes->dm, &levels)); 195c732cbdbSBarry Smith levels++; 196d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 197c732cbdbSBarry Smith } 1989566063dSJacob Faibussowitsch PetscCall(SNESFASSetLevels(snes, levels, NULL)); 19907144faaSPeter Brune fastype = fas->fastype; 2009566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_fas_type", "FAS correction type", "SNESFASSetType", SNESFASTypes, (PetscEnum)fastype, (PetscEnum *)&fastype, &flg)); 2011baa6e33SBarry Smith if (flg) PetscCall(SNESFASSetType(snes, fastype)); 202ee78dd50SPeter Brune 2039566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 2049566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_fas_cycles", "Number of cycles", "SNESFASSetCycles", fas->n_cycles, &m, &flg)); 2051baa6e33SBarry Smith if (flg) PetscCall(SNESFASSetCycles(snes, m)); 2069566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_fas_continuation", "Corrected grid-sequence continuation", "SNESFASSetContinuation", fas->continuation, &continuationflg, &flg)); 2071baa6e33SBarry Smith if (flg) PetscCall(SNESFASSetContinuation(snes, continuationflg)); 208fde0ff24SPeter Brune 2099566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin", "SNESFASSetGalerkin", fas->galerkin, &galerkinflg, &flg)); 2101baa6e33SBarry Smith if (flg) PetscCall(SNESFASSetGalerkin(snes, galerkinflg)); 211ee78dd50SPeter Brune 212928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 2139566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_fas_full_downsweep", "Smooth on the initial down sweep for full FAS cycles", "SNESFASFullSetDownSweep", fas->full_downsweep, &fas->full_downsweep, &flg)); 2149566063dSJacob Faibussowitsch if (flg) PetscCall(SNESFASFullSetDownSweep(snes, fas->full_downsweep)); 2159566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_fas_full_total", "Use total restriction and interpolaton on the indial down and up sweeps for the full FAS cycle", "SNESFASFullSetUseTotal", fas->full_total, &fas->full_total, &flg)); 2169566063dSJacob Faibussowitsch if (flg) PetscCall(SNESFASFullSetTotal(snes, fas->full_total)); 217928e959bSPeter Brune } 218928e959bSPeter Brune 2199566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_fas_smoothup", "Number of post-smoothing steps", "SNESFASSetNumberSmoothUp", fas->max_up_it, &n_up, &upflg)); 220162d76ddSPeter Brune 2219566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_fas_smoothdown", "Number of pre-smoothing steps", "SNESFASSetNumberSmoothDown", fas->max_down_it, &n_down, &downflg)); 222162d76ddSPeter Brune 223d142ab34SLawrence Mitchell { 224d142ab34SLawrence Mitchell PetscViewer viewer; 225d142ab34SLawrence Mitchell PetscViewerFormat format; 2269566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes), ((PetscObject)snes)->options, ((PetscObject)snes)->prefix, "-snes_fas_monitor", &viewer, &format, &monflg)); 227d142ab34SLawrence Mitchell if (monflg) { 228d142ab34SLawrence Mitchell PetscViewerAndFormat *vf; 2299566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer, format, &vf)); 2309566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 2319566063dSJacob Faibussowitsch PetscCall(SNESFASSetMonitor(snes, vf, PETSC_TRUE)); 232d142ab34SLawrence Mitchell } 233d142ab34SLawrence Mitchell } 2340dd27c6cSPeter Brune flg = PETSC_FALSE; 2350dd27c6cSPeter Brune monflg = PETSC_TRUE; 2369566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_fas_log", "Log times for each FAS level", "SNESFASSetLog", monflg, &monflg, &flg)); 2379566063dSJacob Faibussowitsch if (flg) PetscCall(SNESFASSetLog(snes, monflg)); 238ab8d36c9SPeter Brune } 239ee78dd50SPeter Brune 240d0609cedSBarry Smith PetscOptionsHeadEnd(); 241f833ba53SLisandro Dalcin 2428cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 2431baa6e33SBarry Smith if (upflg) PetscCall(SNESFASSetNumberSmoothUp(snes, n_up)); 2441baa6e33SBarry Smith if (downflg) PetscCall(SNESFASSetNumberSmoothDown(snes, n_down)); 245eff52c0eSPeter Brune 2469e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 2479e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 2489e764e56SPeter Brune if (!snes->linesearch) { 2499566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 2509566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2)); 2519e764e56SPeter Brune } 2529e764e56SPeter Brune } 2539e764e56SPeter Brune 254ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 2559566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 2569566063dSJacob Faibussowitsch if (next) PetscCall(SNESSetFromOptions(next)); 257421d9b32SPeter Brune PetscFunctionReturn(0); 258421d9b32SPeter Brune } 259421d9b32SPeter Brune 2609804daf3SBarry Smith #include <petscdraw.h> 2619371c9d4SSatish Balay static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) { 262421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 263656ede7eSPeter Brune PetscBool isFine, iascii, isdraw; 264ab8d36c9SPeter Brune PetscInt i; 265ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 266421d9b32SPeter Brune 267421d9b32SPeter Brune PetscFunctionBegin; 2689566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 269ab8d36c9SPeter Brune if (isFine) { 2709566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 2719566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERDRAW, &isdraw)); 272421d9b32SPeter Brune if (iascii) { 27363a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type is %s, levels=%" PetscInt_FMT ", cycles=%" PetscInt_FMT "\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles)); 274ab8d36c9SPeter Brune if (fas->galerkin) { 2759566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Using Galerkin computed coarse grid function evaluation\n")); 276421d9b32SPeter Brune } else { 2779566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Not using Galerkin computed coarse grid function evaluation\n")); 278421d9b32SPeter Brune } 279ab8d36c9SPeter Brune for (i = 0; i < fas->levels; i++) { 2809566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 2819566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetSmootherUp(levelsnes, &smoothu)); 2829566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetSmootherDown(levelsnes, &smoothd)); 283ab8d36c9SPeter Brune if (!i) { 28463a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n", i)); 285421d9b32SPeter Brune } else { 28663a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 287421d9b32SPeter Brune } 2889566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 289166b3ea4SJed Brown if (smoothd) { 2909566063dSJacob Faibussowitsch PetscCall(SNESView(smoothd, viewer)); 291166b3ea4SJed Brown } else { 2929566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Not yet available\n")); 293166b3ea4SJed Brown } 2949566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 295ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 2969566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Up solver (post-smoother) same as down solver (pre-smoother)\n")); 297ab8d36c9SPeter Brune } else if (i) { 29863a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n", i)); 2999566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 300166b3ea4SJed Brown if (smoothu) { 3019566063dSJacob Faibussowitsch PetscCall(SNESView(smoothu, viewer)); 302166b3ea4SJed Brown } else { 3039566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Not yet available\n")); 304166b3ea4SJed Brown } 3059566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 306ab8d36c9SPeter Brune } 307ab8d36c9SPeter Brune } 308656ede7eSPeter Brune } else if (isdraw) { 309656ede7eSPeter Brune PetscDraw draw; 310b4375e8dSPeter Brune PetscReal x, w, y, bottom, th, wth; 311656ede7eSPeter Brune SNES_FAS *curfas = fas; 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer, 0, &draw)); 3139566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw, &x, &y)); 3149566063dSJacob Faibussowitsch PetscCall(PetscDrawStringGetSize(draw, &wth, &th)); 315656ede7eSPeter Brune bottom = y - th; 316656ede7eSPeter Brune while (curfas) { 317b4375e8dSPeter Brune if (!curfas->smoothu) { 3189566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x, bottom)); 3199566063dSJacob Faibussowitsch if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd, viewer)); 3209566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 321b4375e8dSPeter Brune } else { 322b4375e8dSPeter Brune w = 0.5 * PetscMin(1.0 - x, x); 3239566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x - w, bottom)); 3249566063dSJacob Faibussowitsch if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd, viewer)); 3259566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 3269566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw, x + w, bottom)); 3279566063dSJacob Faibussowitsch if (curfas->smoothu) PetscCall(SNESView(curfas->smoothu, viewer)); 3289566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 329b4375e8dSPeter Brune } 330656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 331656ede7eSPeter Brune bottom -= 5 * th; 3321aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS *)curfas->next->data; 3330298fd71SBarry Smith else curfas = NULL; 334656ede7eSPeter Brune } 335421d9b32SPeter Brune } 336ab8d36c9SPeter Brune } 337421d9b32SPeter Brune PetscFunctionReturn(0); 338421d9b32SPeter Brune } 339421d9b32SPeter Brune 34039bd7f45SPeter Brune /* 34139bd7f45SPeter Brune Defines the action of the downsmoother 34239bd7f45SPeter Brune */ 3439371c9d4SSatish Balay static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) { 344742fe5e2SPeter Brune SNESConvergedReason reason; 345ab8d36c9SPeter Brune Vec FPC; 346ab8d36c9SPeter Brune SNES smoothd; 3476cbb2f26SLawrence Mitchell PetscBool flg; 3480dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 3496e111a19SKarl Rupp 350421d9b32SPeter Brune PetscFunctionBegin; 3519566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetSmootherDown(snes, &smoothd)); 3529566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(smoothd, F)); 3539566063dSJacob Faibussowitsch if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve, smoothd, B, X, 0)); 3549566063dSJacob Faibussowitsch PetscCall(SNESSolve(smoothd, B, X)); 3559566063dSJacob Faibussowitsch if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve, smoothd, B, X, 0)); 356742fe5e2SPeter Brune /* check convergence reason for the smoother */ 3579566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(smoothd, &reason)); 3581ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 359742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 360742fe5e2SPeter Brune PetscFunctionReturn(0); 361742fe5e2SPeter Brune } 3626cbb2f26SLawrence Mitchell 3639566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(smoothd, &FPC, NULL, NULL)); 3649566063dSJacob Faibussowitsch PetscCall(SNESGetAlwaysComputesFinalResidual(smoothd, &flg)); 365*48a46eb9SPierre Jolivet if (!flg) PetscCall(SNESComputeFunction(smoothd, X, FPC)); 3669566063dSJacob Faibussowitsch PetscCall(VecCopy(FPC, F)); 3679566063dSJacob Faibussowitsch if (fnorm) PetscCall(VecNorm(F, NORM_2, fnorm)); 36839bd7f45SPeter Brune PetscFunctionReturn(0); 36939bd7f45SPeter Brune } 37039bd7f45SPeter Brune 37139bd7f45SPeter Brune /* 37207144faaSPeter Brune Defines the action of the upsmoother 37339bd7f45SPeter Brune */ 3749371c9d4SSatish Balay static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) { 37539bd7f45SPeter Brune SNESConvergedReason reason; 376ab8d36c9SPeter Brune Vec FPC; 377ab8d36c9SPeter Brune SNES smoothu; 3786cbb2f26SLawrence Mitchell PetscBool flg; 3790dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 380ab8d36c9SPeter Brune 3816e111a19SKarl Rupp PetscFunctionBegin; 3829566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetSmootherUp(snes, &smoothu)); 3839566063dSJacob Faibussowitsch if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve, smoothu, 0, 0, 0)); 3849566063dSJacob Faibussowitsch PetscCall(SNESSolve(smoothu, B, X)); 3859566063dSJacob Faibussowitsch if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve, smoothu, 0, 0, 0)); 38639bd7f45SPeter Brune /* check convergence reason for the smoother */ 3879566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(smoothu, &reason)); 3881ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 38939bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 39039bd7f45SPeter Brune PetscFunctionReturn(0); 39139bd7f45SPeter Brune } 3929566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(smoothu, &FPC, NULL, NULL)); 3939566063dSJacob Faibussowitsch PetscCall(SNESGetAlwaysComputesFinalResidual(smoothu, &flg)); 394*48a46eb9SPierre Jolivet if (!flg) PetscCall(SNESComputeFunction(smoothu, X, FPC)); 3959566063dSJacob Faibussowitsch PetscCall(VecCopy(FPC, F)); 3969566063dSJacob Faibussowitsch if (fnorm) PetscCall(VecNorm(F, NORM_2, fnorm)); 39739bd7f45SPeter Brune PetscFunctionReturn(0); 39839bd7f45SPeter Brune } 39939bd7f45SPeter Brune 400938e4a01SJed Brown /*@ 401938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 402938e4a01SJed Brown 403938e4a01SJed Brown Collective 404938e4a01SJed Brown 4054165533cSJose E. Roman Input Parameter: 406938e4a01SJed Brown . snes - SNESFAS 407938e4a01SJed Brown 4084165533cSJose E. Roman Output Parameter: 409938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 410938e4a01SJed Brown 411938e4a01SJed Brown Level: developer 412938e4a01SJed Brown 413db781477SPatrick Sanan .seealso: `SNESFASSetRestriction()`, `SNESFASRestrict()` 414938e4a01SJed Brown @*/ 4159371c9d4SSatish Balay PetscErrorCode SNESFASCreateCoarseVec(SNES snes, Vec *Xcoarse) { 416f833ba53SLisandro Dalcin SNES_FAS *fas; 417938e4a01SJed Brown 418938e4a01SJed Brown PetscFunctionBegin; 419f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS); 420f833ba53SLisandro Dalcin PetscValidPointer(Xcoarse, 2); 421f833ba53SLisandro Dalcin fas = (SNES_FAS *)snes->data; 4221aa26658SKarl Rupp if (fas->rscale) { 4239566063dSJacob Faibussowitsch PetscCall(VecDuplicate(fas->rscale, Xcoarse)); 424f5af7f23SKarl Rupp } else if (fas->inject) { 4259566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(fas->inject, Xcoarse, NULL)); 42613903a91SSatish Balay } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Must set restriction or injection"); 427938e4a01SJed Brown PetscFunctionReturn(0); 428938e4a01SJed Brown } 429938e4a01SJed Brown 430e9923e8dSJed Brown /*@ 431e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 432e9923e8dSJed Brown 433e9923e8dSJed Brown Collective 434e9923e8dSJed Brown 4354165533cSJose E. Roman Input Parameters: 436e9923e8dSJed Brown + fine - SNES from which to restrict 437e9923e8dSJed Brown - Xfine - vector to restrict 438e9923e8dSJed Brown 4394165533cSJose E. Roman Output Parameter: 440e9923e8dSJed Brown . Xcoarse - result of restriction 441e9923e8dSJed Brown 442e9923e8dSJed Brown Level: developer 443e9923e8dSJed Brown 444db781477SPatrick Sanan .seealso: `SNESFASSetRestriction()`, `SNESFASSetInjection()` 445e9923e8dSJed Brown @*/ 4469371c9d4SSatish Balay PetscErrorCode SNESFASRestrict(SNES fine, Vec Xfine, Vec Xcoarse) { 447f833ba53SLisandro Dalcin SNES_FAS *fas; 448e9923e8dSJed Brown 449e9923e8dSJed Brown PetscFunctionBegin; 450f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(fine, SNES_CLASSID, 1, SNESFAS); 451e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine, VEC_CLASSID, 2); 452e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse, VEC_CLASSID, 3); 453f833ba53SLisandro Dalcin fas = (SNES_FAS *)fine->data; 454e9923e8dSJed Brown if (fas->inject) { 4559566063dSJacob Faibussowitsch PetscCall(MatRestrict(fas->inject, Xfine, Xcoarse)); 456e9923e8dSJed Brown } else { 4579566063dSJacob Faibussowitsch PetscCall(MatRestrict(fas->restrct, Xfine, Xcoarse)); 4589566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Xcoarse, fas->rscale, Xcoarse)); 459e9923e8dSJed Brown } 460e9923e8dSJed Brown PetscFunctionReturn(0); 461e9923e8dSJed Brown } 462e9923e8dSJed Brown 46339bd7f45SPeter Brune /* 46439bd7f45SPeter Brune 4656dfbafefSToby Isaac Performs a variant of FAS using the interpolated total coarse solution 4666dfbafefSToby Isaac 4676dfbafefSToby Isaac fine problem: F(x) = b 4686dfbafefSToby Isaac coarse problem: F^c(x^c) = Rb, Initial guess Rx 4696dfbafefSToby Isaac interpolated solution: x^f = I x^c (total solution interpolation 4706dfbafefSToby Isaac 4716dfbafefSToby Isaac */ 4729371c9d4SSatish Balay static PetscErrorCode SNESFASInterpolatedCoarseSolution(SNES snes, Vec X, Vec X_new) { 4736dfbafefSToby Isaac Vec X_c, B_c; 4746dfbafefSToby Isaac SNESConvergedReason reason; 4756dfbafefSToby Isaac SNES next; 4766dfbafefSToby Isaac Mat restrct, interpolate; 4776dfbafefSToby Isaac SNES_FAS *fasc; 4786dfbafefSToby Isaac 4796dfbafefSToby Isaac PetscFunctionBegin; 4809566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 4816dfbafefSToby Isaac if (next) { 4826dfbafefSToby Isaac fasc = (SNES_FAS *)next->data; 4836dfbafefSToby Isaac 4849566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetRestriction(snes, &restrct)); 4859566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate)); 4866dfbafefSToby Isaac 4876dfbafefSToby Isaac X_c = next->vec_sol; 4886dfbafefSToby Isaac 4899566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0)); 4906dfbafefSToby Isaac /* restrict the total solution: Rb */ 4919566063dSJacob Faibussowitsch PetscCall(SNESFASRestrict(snes, X, X_c)); 4926dfbafefSToby Isaac B_c = next->vec_rhs; 4936dfbafefSToby Isaac if (snes->vec_rhs) { 4946dfbafefSToby Isaac /* restrict the total rhs defect: Rb */ 4959566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, snes->vec_rhs, B_c)); 4966dfbafefSToby Isaac } else { 4979566063dSJacob Faibussowitsch PetscCall(VecSet(B_c, 0.)); 4986dfbafefSToby Isaac } 4999566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0)); 5006dfbafefSToby Isaac 5019566063dSJacob Faibussowitsch PetscCall(SNESSolve(next, B_c, X_c)); 5029566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next, &reason)); 5036dfbafefSToby Isaac if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 5046dfbafefSToby Isaac snes->reason = SNES_DIVERGED_INNER; 5056dfbafefSToby Isaac PetscFunctionReturn(0); 5066dfbafefSToby Isaac } 5076dfbafefSToby Isaac /* x^f <- Ix^c*/ 5086dfbafefSToby Isaac DM dmc, dmf; 5096dfbafefSToby Isaac 5109566063dSJacob Faibussowitsch PetscCall(SNESGetDM(next, &dmc)); 5119566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dmf)); 5129566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0)); 5139566063dSJacob Faibussowitsch PetscCall(DMInterpolateSolution(dmc, dmf, interpolate, X_c, X_new)); 5149566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0)); 5159566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)X_c, "Coarse solution")); 5169566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view")); 5179566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)X_new, "Updated Fine solution")); 5189566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view")); 5196dfbafefSToby Isaac } 5206dfbafefSToby Isaac PetscFunctionReturn(0); 5216dfbafefSToby Isaac } 5226dfbafefSToby Isaac 5236dfbafefSToby Isaac /* 5246dfbafefSToby Isaac 52539bd7f45SPeter Brune Performs the FAS coarse correction as: 52639bd7f45SPeter Brune 527b5527d98SMatthew G. Knepley fine problem: F(x) = b 528b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c 52939bd7f45SPeter Brune 530b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 53139bd7f45SPeter Brune 53239bd7f45SPeter Brune */ 5339371c9d4SSatish Balay PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) { 53439bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 53539bd7f45SPeter Brune SNESConvergedReason reason; 536ab8d36c9SPeter Brune SNES next; 537ab8d36c9SPeter Brune Mat restrct, interpolate; 5380dd27c6cSPeter Brune SNES_FAS *fasc; 5395fd66863SKarl Rupp 54039bd7f45SPeter Brune PetscFunctionBegin; 5419566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 542ab8d36c9SPeter Brune if (next) { 5430dd27c6cSPeter Brune fasc = (SNES_FAS *)next->data; 5440dd27c6cSPeter Brune 5459566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetRestriction(snes, &restrct)); 5469566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate)); 547ab8d36c9SPeter Brune 548ab8d36c9SPeter Brune X_c = next->vec_sol; 549ab8d36c9SPeter Brune Xo_c = next->work[0]; 550ab8d36c9SPeter Brune F_c = next->vec_func; 551ab8d36c9SPeter Brune B_c = next->vec_rhs; 552efe1f98aSPeter Brune 5539566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0)); 5549566063dSJacob Faibussowitsch PetscCall(SNESFASRestrict(snes, X, Xo_c)); 5555609cbf2SMatthew G. Knepley /* restrict the defect: R(F(x) - b) */ 5569566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, F, B_c)); 5579566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0)); 5580dd27c6cSPeter Brune 5599566063dSJacob Faibussowitsch if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual, next, 0, 0, 0)); 5605609cbf2SMatthew G. Knepley /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */ 5619566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(next, Xo_c, F_c)); 5629566063dSJacob Faibussowitsch if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual, next, 0, 0, 0)); 5630dd27c6cSPeter Brune 5640dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 5659566063dSJacob Faibussowitsch PetscCall(VecCopy(B_c, X_c)); 5669566063dSJacob Faibussowitsch PetscCall(VecCopy(F_c, B_c)); 5679566063dSJacob Faibussowitsch PetscCall(VecCopy(X_c, F_c)); 568ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 5699566063dSJacob Faibussowitsch PetscCall(VecCopy(Xo_c, X_c)); 570c90fad12SPeter Brune 571c90fad12SPeter Brune /* recurse to the next level */ 5729566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next, F_c)); 5739566063dSJacob Faibussowitsch PetscCall(SNESSolve(next, B_c, X_c)); 5749566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next, &reason)); 575742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 576742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 577742fe5e2SPeter Brune PetscFunctionReturn(0); 578742fe5e2SPeter Brune } 579fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 5809566063dSJacob Faibussowitsch PetscCall(VecAXPY(X_c, -1.0, Xo_c)); 5810dd27c6cSPeter Brune 5829566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict, snes, 0, 0, 0)); 5839566063dSJacob Faibussowitsch PetscCall(MatInterpolateAdd(interpolate, X_c, X, X_new)); 5849566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict, snes, 0, 0, 0)); 5859566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)X_c, "Coarse correction")); 5869566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view")); 5879566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject)X_new, "Updated Fine solution")); 5889566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view")); 589293a7e31SPeter Brune } 59039bd7f45SPeter Brune PetscFunctionReturn(0); 59139bd7f45SPeter Brune } 59239bd7f45SPeter Brune 59339bd7f45SPeter Brune /* 59439bd7f45SPeter Brune 59539bd7f45SPeter Brune The additive cycle looks like: 59639bd7f45SPeter Brune 59707144faaSPeter Brune xhat = x 59807144faaSPeter Brune xhat = dS(x, b) 59907144faaSPeter Brune x = coarsecorrection(xhat, b_d) 60007144faaSPeter Brune x = x + nu*(xhat - x); 60139bd7f45SPeter Brune (optional) x = uS(x, b) 60239bd7f45SPeter Brune 60339bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 60439bd7f45SPeter Brune 60539bd7f45SPeter Brune */ 6069371c9d4SSatish Balay static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) { 60707144faaSPeter Brune Vec F, B, Xhat; 60822c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 60907144faaSPeter Brune SNESConvergedReason reason; 61022c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 611422a814eSBarry Smith SNESLineSearchReason lsresult; 612ab8d36c9SPeter Brune SNES next; 613ab8d36c9SPeter Brune Mat restrct, interpolate; 6140dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data, *fasc; 6150dd27c6cSPeter Brune 61639bd7f45SPeter Brune PetscFunctionBegin; 6179566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 61839bd7f45SPeter Brune F = snes->vec_func; 61939bd7f45SPeter Brune B = snes->vec_rhs; 620e7f468e7SPeter Brune Xhat = snes->work[1]; 6219566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xhat)); 62207144faaSPeter Brune /* recurse first */ 623ab8d36c9SPeter Brune if (next) { 6240dd27c6cSPeter Brune fasc = (SNES_FAS *)next->data; 6259566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetRestriction(snes, &restrct)); 6269566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate)); 6279566063dSJacob Faibussowitsch if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual, snes, 0, 0, 0)); 6289566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, Xhat, F)); 6299566063dSJacob Faibussowitsch if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual, snes, 0, 0, 0)); 6309566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 631ab8d36c9SPeter Brune X_c = next->vec_sol; 632ab8d36c9SPeter Brune Xo_c = next->work[0]; 633ab8d36c9SPeter Brune F_c = next->vec_func; 634ab8d36c9SPeter Brune B_c = next->vec_rhs; 63539bd7f45SPeter Brune 6369566063dSJacob Faibussowitsch PetscCall(SNESFASRestrict(snes, Xhat, Xo_c)); 63707144faaSPeter Brune /* restrict the defect */ 6389566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, F, B_c)); 63907144faaSPeter Brune 64007144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 6419566063dSJacob Faibussowitsch if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual, next, 0, 0, 0)); 6429566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(next, Xo_c, F_c)); 6439566063dSJacob Faibussowitsch if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual, next, 0, 0, 0)); 6449566063dSJacob Faibussowitsch PetscCall(VecCopy(B_c, X_c)); 6459566063dSJacob Faibussowitsch PetscCall(VecCopy(F_c, B_c)); 6469566063dSJacob Faibussowitsch PetscCall(VecCopy(X_c, F_c)); 64707144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 6489566063dSJacob Faibussowitsch PetscCall(VecCopy(Xo_c, X_c)); 64907144faaSPeter Brune 65007144faaSPeter Brune /* recurse */ 6519566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next, F_c)); 6529566063dSJacob Faibussowitsch PetscCall(SNESSolve(next, B_c, X_c)); 65307144faaSPeter Brune 65407144faaSPeter Brune /* smooth on this level */ 6559566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &fnorm)); 65607144faaSPeter Brune 6579566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next, &reason)); 65807144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 65907144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 66007144faaSPeter Brune PetscFunctionReturn(0); 66107144faaSPeter Brune } 66207144faaSPeter Brune 66307144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 6649566063dSJacob Faibussowitsch PetscCall(VecAYPX(X_c, -1.0, Xo_c)); 6659566063dSJacob Faibussowitsch PetscCall(MatInterpolate(interpolate, X_c, Xhat)); 66607144faaSPeter Brune 667ddebd997SPeter Brune /* additive correction of the coarse direction*/ 6689566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat)); 6699566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult)); 6709566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm)); 671422a814eSBarry Smith if (lsresult) { 6729e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 6739e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 6749e764e56SPeter Brune PetscFunctionReturn(0); 6759e764e56SPeter Brune } 6769e764e56SPeter Brune } 67707144faaSPeter Brune } else { 6789566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 67907144faaSPeter Brune } 68039bd7f45SPeter Brune PetscFunctionReturn(0); 68139bd7f45SPeter Brune } 68239bd7f45SPeter Brune 68339bd7f45SPeter Brune /* 68439bd7f45SPeter Brune 68539bd7f45SPeter Brune Defines the FAS cycle as: 68639bd7f45SPeter Brune 687b5527d98SMatthew G. Knepley fine problem: F(x) = b 68839bd7f45SPeter Brune coarse problem: F^c(x) = b^c 68939bd7f45SPeter Brune 690b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 69139bd7f45SPeter Brune 69239bd7f45SPeter Brune correction: 69339bd7f45SPeter Brune 69439bd7f45SPeter Brune x = x + I(x^c - Rx) 69539bd7f45SPeter Brune 69639bd7f45SPeter Brune */ 6979371c9d4SSatish Balay static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) { 69839bd7f45SPeter Brune Vec F, B; 69934d65b3cSPeter Brune SNES next; 70039bd7f45SPeter Brune 70139bd7f45SPeter Brune PetscFunctionBegin; 70239bd7f45SPeter Brune F = snes->vec_func; 70339bd7f45SPeter Brune B = snes->vec_rhs; 70439bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 7059566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 7069566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 70734d65b3cSPeter Brune if (next) { 7089566063dSJacob Faibussowitsch PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 7099566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 710fe6f9142SPeter Brune } 711fa9694d7SPeter Brune PetscFunctionReturn(0); 712421d9b32SPeter Brune } 713421d9b32SPeter Brune 7149371c9d4SSatish Balay static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) { 7158c94862eSPeter Brune SNES next; 7168c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 7178c94862eSPeter Brune PetscBool isFine; 7188c94862eSPeter Brune 7198c94862eSPeter Brune PetscFunctionBegin; 7208c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 7219566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 7229566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 7238c94862eSPeter Brune fas->full_stage = 0; 7249566063dSJacob Faibussowitsch if (next) PetscCall(SNESFASCycleSetupPhase_Full(next)); 7258c94862eSPeter Brune PetscFunctionReturn(0); 7268c94862eSPeter Brune } 7278c94862eSPeter Brune 7289371c9d4SSatish Balay static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) { 7298c94862eSPeter Brune Vec F, B; 7308c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS *)snes->data; 7318c94862eSPeter Brune PetscBool isFine; 7328c94862eSPeter Brune SNES next; 7338c94862eSPeter Brune 7348c94862eSPeter Brune PetscFunctionBegin; 7358c94862eSPeter Brune F = snes->vec_func; 7368c94862eSPeter Brune B = snes->vec_rhs; 7379566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 7389566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 7398c94862eSPeter Brune 7401baa6e33SBarry Smith if (isFine) PetscCall(SNESFASCycleSetupPhase_Full(snes)); 7418c94862eSPeter Brune 7428c94862eSPeter Brune if (fas->full_stage == 0) { 743928e959bSPeter Brune /* downsweep */ 7448c94862eSPeter Brune if (next) { 7458c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 7469566063dSJacob Faibussowitsch if (fas->full_downsweep) PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 747e140b65aSPatrick Farrell fas->full_downsweep = PETSC_TRUE; 7489566063dSJacob Faibussowitsch if (fas->full_total) PetscCall(SNESFASInterpolatedCoarseSolution(snes, X, X)); 7499566063dSJacob Faibussowitsch else PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 7506dfbafefSToby Isaac fas->full_total = PETSC_FALSE; 7519566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 7528c94862eSPeter Brune if (fas->level != 1) next->max_its -= 1; 7538c94862eSPeter Brune } else { 754a3a80b83SMatthew G. Knepley /* The smoother on the coarse level is the coarse solver */ 7559566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 7568c94862eSPeter Brune } 7578c94862eSPeter Brune fas->full_stage = 1; 7588c94862eSPeter Brune } else if (fas->full_stage == 1) { 7599566063dSJacob Faibussowitsch if (snes->iter == 0) PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 7608c94862eSPeter Brune if (next) { 7619566063dSJacob Faibussowitsch PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 7629566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 7638c94862eSPeter Brune } 7648c94862eSPeter Brune } 7658c94862eSPeter Brune /* final v-cycle */ 7668c94862eSPeter Brune if (isFine) { 7678c94862eSPeter Brune if (next) { 7689566063dSJacob Faibussowitsch PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 7699566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 7708c94862eSPeter Brune } 7718c94862eSPeter Brune } 7728c94862eSPeter Brune PetscFunctionReturn(0); 7738c94862eSPeter Brune } 7748c94862eSPeter Brune 7759371c9d4SSatish Balay static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) { 77634d65b3cSPeter Brune Vec F, B; 77734d65b3cSPeter Brune SNES next; 77834d65b3cSPeter Brune 77934d65b3cSPeter Brune PetscFunctionBegin; 78034d65b3cSPeter Brune F = snes->vec_func; 78134d65b3cSPeter Brune B = snes->vec_rhs; 7829566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 78334d65b3cSPeter Brune if (next) { 7849566063dSJacob Faibussowitsch PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 7859566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 78634d65b3cSPeter Brune } else { 7879566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 78834d65b3cSPeter Brune } 78934d65b3cSPeter Brune PetscFunctionReturn(0); 79034d65b3cSPeter Brune } 79134d65b3cSPeter Brune 792fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE; 793fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n" 794fffbeea8SBarry Smith " title = {Composing Scalable Nonlinear Algebraic Solvers},\n" 795fffbeea8SBarry Smith " author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n" 796fffbeea8SBarry Smith " year = 2013,\n" 797fffbeea8SBarry Smith " type = Preprint,\n" 798fffbeea8SBarry Smith " number = {ANL/MCS-P2010-0112},\n" 799fffbeea8SBarry Smith " institution = {Argonne National Laboratory}\n}\n"; 800fffbeea8SBarry Smith 8019371c9d4SSatish Balay static PetscErrorCode SNESSolve_FAS(SNES snes) { 802f833ba53SLisandro Dalcin PetscInt i; 803ddb5aff1SPeter Brune Vec X, F; 804fe6f9142SPeter Brune PetscReal fnorm; 805b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS *)snes->data, *ffas; 806b17ce1afSJed Brown DM dm; 807e70c42e5SPeter Brune PetscBool isFine; 808b17ce1afSJed Brown 809421d9b32SPeter Brune PetscFunctionBegin; 8100b121fc5SBarry Smith PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 811c579b300SPatrick Farrell 8129566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite)); 813fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 814fa9694d7SPeter Brune X = snes->vec_sol; 815f5a6d4f9SBarry Smith F = snes->vec_func; 816293a7e31SPeter Brune 8179566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 818293a7e31SPeter Brune /* norm setup */ 8199566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 820fe6f9142SPeter Brune snes->iter = 0; 821f833ba53SLisandro Dalcin snes->norm = 0; 8229566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 823e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 8249566063dSJacob Faibussowitsch if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual, snes, 0, 0, 0)); 8259566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, X, F)); 8269566063dSJacob Faibussowitsch if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual, snes, 0, 0, 0)); 8271aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 828e4ed7901SPeter Brune 8299566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 830422a814eSBarry Smith SNESCheckFunctionNorm(snes, fnorm); 8319566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 832fe6f9142SPeter Brune snes->norm = fnorm; 8339566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 8349566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, fnorm, 0)); 8359566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, fnorm)); 836fe6f9142SPeter Brune 837fe6f9142SPeter Brune /* test convergence */ 838dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, 0, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP); 839fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 840e4ed7901SPeter Brune 841b9c2fdf1SPeter Brune if (isFine) { 842b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 8439566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dm)); 844b17ce1afSJed Brown for (ffas = fas; ffas->next; ffas = (SNES_FAS *)ffas->next->data) { 845b17ce1afSJed Brown DM dmcoarse; 8469566063dSJacob Faibussowitsch PetscCall(SNESGetDM(ffas->next, &dmcoarse)); 8479566063dSJacob Faibussowitsch PetscCall(DMRestrict(dm, ffas->restrct, ffas->rscale, ffas->inject, dmcoarse)); 848b17ce1afSJed Brown dm = dmcoarse; 849b17ce1afSJed Brown } 850b9c2fdf1SPeter Brune } 851b17ce1afSJed Brown 852f833ba53SLisandro Dalcin for (i = 0; i < snes->max_its; i++) { 853fe6f9142SPeter Brune /* Call general purpose update function */ 854dbbe0bcdSBarry Smith PetscTryTypeMethod(snes, update, snes->iter); 855f833ba53SLisandro Dalcin 85607144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 8579566063dSJacob Faibussowitsch PetscCall(SNESFASCycle_Multiplicative(snes, X)); 8588c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_ADDITIVE) { 8599566063dSJacob Faibussowitsch PetscCall(SNESFASCycle_Additive(snes, X)); 8608c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_FULL) { 8619566063dSJacob Faibussowitsch PetscCall(SNESFASCycle_Full(snes, X)); 86234d65b3cSPeter Brune } else if (fas->fastype == SNES_FAS_KASKADE) { 8639566063dSJacob Faibussowitsch PetscCall(SNESFASCycle_Kaskade(snes, X)); 8646c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "Unsupported FAS type"); 865742fe5e2SPeter Brune 866742fe5e2SPeter Brune /* check for FAS cycle divergence */ 8671aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 868b9c2fdf1SPeter Brune 869c90fad12SPeter Brune /* Monitor convergence */ 8709566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 871c90fad12SPeter Brune snes->iter = i + 1; 8729566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 8739566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0)); 8749566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes, snes->iter, snes->norm)); 875c90fad12SPeter Brune /* Test for convergence */ 87666585501SPeter Brune if (isFine) { 877dbbe0bcdSBarry Smith PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, snes->norm, &snes->reason, snes->cnvP); 878c90fad12SPeter Brune if (snes->reason) break; 879fe6f9142SPeter Brune } 88066585501SPeter Brune } 881f833ba53SLisandro Dalcin if (i == snes->max_its) { 88263a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", i)); 883fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 884fe6f9142SPeter Brune } 885421d9b32SPeter Brune PetscFunctionReturn(0); 886421d9b32SPeter Brune } 88740244768SBarry Smith 88840244768SBarry Smith /*MC 88940244768SBarry Smith 89040244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 89140244768SBarry Smith 89240244768SBarry Smith The nonlinear problem is solved by correction using coarse versions 89340244768SBarry Smith of the nonlinear problem. This problem is perturbed so that a projected 89440244768SBarry Smith solution of the fine problem elicits no correction from the coarse problem. 89540244768SBarry Smith 89640244768SBarry Smith Options Database: 89740244768SBarry Smith + -snes_fas_levels - The number of levels 89840244768SBarry Smith . -snes_fas_cycles<1> - The number of cycles -- 1 for V, 2 for W 89940244768SBarry Smith . -snes_fas_type<additive,multiplicative,full,kaskade> - Additive or multiplicative cycle 90040244768SBarry Smith . -snes_fas_galerkin<PETSC_FALSE> - Form coarse problems by projection back upon the fine problem 90140244768SBarry Smith . -snes_fas_smoothup<1> - The number of iterations of the post-smoother 90240244768SBarry Smith . -snes_fas_smoothdown<1> - The number of iterations of the pre-smoother 90340244768SBarry Smith . -snes_fas_monitor - Monitor progress of all of the levels 90440244768SBarry Smith . -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS 90540244768SBarry Smith . -fas_levels_snes_ - SNES options for all smoothers 90640244768SBarry Smith . -fas_levels_cycle_snes_ - SNES options for all cycles 90740244768SBarry Smith . -fas_levels_i_snes_ - SNES options for the smoothers on level i 90840244768SBarry Smith . -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i 90940244768SBarry Smith - -fas_coarse_snes_ - SNES options for the coarsest smoother 91040244768SBarry Smith 91140244768SBarry Smith Notes: 91240244768SBarry Smith The organization of the FAS solver is slightly different from the organization of PCMG 91340244768SBarry Smith As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance. 91440244768SBarry Smith The cycle SNES instance may be used for monitoring convergence on a particular level. 91540244768SBarry Smith 91640244768SBarry Smith Level: beginner 91740244768SBarry Smith 91840244768SBarry Smith References: 919606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 92040244768SBarry Smith SIAM Review, 57(4), 2015 92140244768SBarry Smith 922db781477SPatrick Sanan .seealso: `PCMG`, `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESType` 92340244768SBarry Smith M*/ 92440244768SBarry Smith 9259371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes) { 92640244768SBarry Smith SNES_FAS *fas; 92740244768SBarry Smith 92840244768SBarry Smith PetscFunctionBegin; 92940244768SBarry Smith snes->ops->destroy = SNESDestroy_FAS; 93040244768SBarry Smith snes->ops->setup = SNESSetUp_FAS; 93140244768SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_FAS; 93240244768SBarry Smith snes->ops->view = SNESView_FAS; 93340244768SBarry Smith snes->ops->solve = SNESSolve_FAS; 93440244768SBarry Smith snes->ops->reset = SNESReset_FAS; 93540244768SBarry Smith 93640244768SBarry Smith snes->usesksp = PETSC_FALSE; 937efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 93840244768SBarry Smith 93940244768SBarry Smith if (!snes->tolerancesset) { 94040244768SBarry Smith snes->max_funcs = 30000; 94140244768SBarry Smith snes->max_its = 10000; 94240244768SBarry Smith } 94340244768SBarry Smith 9444fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 9454fc747eaSLawrence Mitchell 9469566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes, &fas)); 94740244768SBarry Smith 94840244768SBarry Smith snes->data = (void *)fas; 94940244768SBarry Smith fas->level = 0; 95040244768SBarry Smith fas->levels = 1; 95140244768SBarry Smith fas->n_cycles = 1; 95240244768SBarry Smith fas->max_up_it = 1; 95340244768SBarry Smith fas->max_down_it = 1; 95440244768SBarry Smith fas->smoothu = NULL; 95540244768SBarry Smith fas->smoothd = NULL; 95640244768SBarry Smith fas->next = NULL; 95740244768SBarry Smith fas->previous = NULL; 95840244768SBarry Smith fas->fine = snes; 95940244768SBarry Smith fas->interpolate = NULL; 96040244768SBarry Smith fas->restrct = NULL; 96140244768SBarry Smith fas->inject = NULL; 96240244768SBarry Smith fas->usedmfornumberoflevels = PETSC_FALSE; 96340244768SBarry Smith fas->fastype = SNES_FAS_MULTIPLICATIVE; 96440244768SBarry Smith fas->full_downsweep = PETSC_FALSE; 9656dfbafefSToby Isaac fas->full_total = PETSC_FALSE; 96640244768SBarry Smith 96740244768SBarry Smith fas->eventsmoothsetup = 0; 96840244768SBarry Smith fas->eventsmoothsolve = 0; 96940244768SBarry Smith fas->eventresidual = 0; 97040244768SBarry Smith fas->eventinterprestrict = 0; 97140244768SBarry Smith PetscFunctionReturn(0); 97240244768SBarry Smith } 973