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 640244768SBarry Smith static PetscErrorCode SNESReset_FAS(SNES snes) 7421d9b32SPeter Brune { 8421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 9421d9b32SPeter Brune 10421d9b32SPeter Brune PetscFunctionBegin; 119566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&fas->smoothu)); 129566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&fas->smoothd)); 139566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->inject)); 149566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->interpolate)); 159566063dSJacob Faibussowitsch PetscCall(MatDestroy(&fas->restrct)); 169566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->rscale)); 179566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->Xg)); 189566063dSJacob Faibussowitsch PetscCall(VecDestroy(&fas->Fg)); 199566063dSJacob Faibussowitsch if (fas->next) PetscCall(SNESReset(fas->next)); 20421d9b32SPeter Brune PetscFunctionReturn(0); 21421d9b32SPeter Brune } 22421d9b32SPeter Brune 2340244768SBarry Smith static PetscErrorCode SNESDestroy_FAS(SNES snes) 24421d9b32SPeter Brune { 25421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 26421d9b32SPeter Brune 27421d9b32SPeter Brune PetscFunctionBegin; 28421d9b32SPeter Brune /* recursively resets and then destroys */ 299566063dSJacob Faibussowitsch PetscCall(SNESReset_FAS(snes)); 309566063dSJacob Faibussowitsch PetscCall(SNESDestroy(&fas->next)); 319566063dSJacob Faibussowitsch PetscCall(PetscFree(fas)); 32421d9b32SPeter Brune PetscFunctionReturn(0); 33421d9b32SPeter Brune } 34421d9b32SPeter Brune 35f833ba53SLisandro Dalcin static PetscErrorCode SNESFASSetUpLineSearch_Private(SNES snes, SNES smooth) 36f833ba53SLisandro Dalcin { 37f833ba53SLisandro Dalcin SNESLineSearch linesearch; 38f833ba53SLisandro Dalcin SNESLineSearch slinesearch; 39f833ba53SLisandro Dalcin void *lsprectx,*lspostctx; 40f833ba53SLisandro Dalcin PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 41f833ba53SLisandro Dalcin PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 42f833ba53SLisandro Dalcin 43f833ba53SLisandro Dalcin PetscFunctionBegin; 44f833ba53SLisandro Dalcin if (!snes->linesearch) PetscFunctionReturn(0); 459566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes,&linesearch)); 469566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(smooth,&slinesearch)); 479566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx)); 489566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx)); 499566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx)); 509566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx)); 519566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch)); 52f833ba53SLisandro Dalcin PetscFunctionReturn(0); 53f833ba53SLisandro Dalcin } 54f833ba53SLisandro Dalcin 55f833ba53SLisandro Dalcin static PetscErrorCode SNESFASCycleSetUpSmoother_Private(SNES snes, SNES smooth) 56f833ba53SLisandro Dalcin { 57f833ba53SLisandro Dalcin SNES_FAS *fas = (SNES_FAS*) snes->data; 58f833ba53SLisandro Dalcin 59f833ba53SLisandro Dalcin PetscFunctionBegin; 609566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)smooth)); 619566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(smooth)); 629566063dSJacob Faibussowitsch PetscCall(SNESFASSetUpLineSearch_Private(snes, smooth)); 63f833ba53SLisandro Dalcin 649566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)snes->vec_sol)); 659566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)snes->vec_sol_update)); 669566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)snes->vec_func)); 67f833ba53SLisandro Dalcin smooth->vec_sol = snes->vec_sol; 68f833ba53SLisandro Dalcin smooth->vec_sol_update = snes->vec_sol_update; 69f833ba53SLisandro Dalcin smooth->vec_func = snes->vec_func; 70f833ba53SLisandro Dalcin 719566063dSJacob Faibussowitsch if (fas->eventsmoothsetup) PetscCall(PetscLogEventBegin(fas->eventsmoothsetup,smooth,0,0,0)); 729566063dSJacob Faibussowitsch PetscCall(SNESSetUp(smooth)); 739566063dSJacob Faibussowitsch if (fas->eventsmoothsetup) PetscCall(PetscLogEventEnd(fas->eventsmoothsetup,smooth,0,0,0)); 74f833ba53SLisandro Dalcin PetscFunctionReturn(0); 75f833ba53SLisandro Dalcin } 76f833ba53SLisandro Dalcin 7740244768SBarry Smith static PetscErrorCode SNESSetUp_FAS(SNES snes) 78421d9b32SPeter Brune { 7948bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 80d1adcc6fSPeter Brune PetscInt dm_levels; 81ab8d36c9SPeter Brune SNES next; 82f833ba53SLisandro Dalcin PetscBool isFine, hasCreateRestriction, hasCreateInjection; 83eff52c0eSPeter Brune 846b2b7091SBarry Smith PetscFunctionBegin; 859566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 86ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 879566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(snes->dm,&dm_levels)); 88d1adcc6fSPeter Brune dm_levels++; 89cc05f883SPeter Brune if (dm_levels > fas->levels) { 903dccd265SPeter Brune /* reset the number of levels */ 919566063dSJacob Faibussowitsch PetscCall(SNESFASSetLevels(snes,dm_levels,NULL)); 929566063dSJacob Faibussowitsch PetscCall(SNESSetFromOptions(snes)); 93d1adcc6fSPeter Brune } 94d1adcc6fSPeter Brune } 959566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 96ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 973dccd265SPeter Brune 989566063dSJacob Faibussowitsch PetscCall(SNESSetWorkVecs(snes, 2)); /* work vectors used for intergrid transfers */ 99cc05f883SPeter Brune 100ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 101ab8d36c9SPeter Brune if (!fas->smoothd) { 1029566063dSJacob Faibussowitsch PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd)); 103ab8d36c9SPeter Brune } 104ab8d36c9SPeter Brune 10579d9a41aSPeter Brune if (snes->dm) { 106ab8d36c9SPeter Brune /* set the smoother DMs properly */ 1079566063dSJacob Faibussowitsch if (fas->smoothu) PetscCall(SNESSetDM(fas->smoothu, snes->dm)); 1089566063dSJacob Faibussowitsch PetscCall(SNESSetDM(fas->smoothd, snes->dm)); 10979d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 110ab8d36c9SPeter Brune if (next) { 11179d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 112ab8d36c9SPeter Brune if (!next->dm) { 1139566063dSJacob Faibussowitsch PetscCall(DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm)); 1149566063dSJacob Faibussowitsch PetscCall(SNESSetDM(next, next->dm)); 11579d9a41aSPeter Brune } 11679d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 11779d9a41aSPeter Brune if (!fas->interpolate) { 1189566063dSJacob Faibussowitsch PetscCall(DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale)); 119bccf9bb3SJed Brown if (!fas->restrct) { 1209566063dSJacob Faibussowitsch PetscCall(DMHasCreateRestriction(next->dm, &hasCreateRestriction)); 1210a7266b2SPatrick Farrell /* DM can create restrictions, use that */ 1220a7266b2SPatrick Farrell if (hasCreateRestriction) { 1239566063dSJacob Faibussowitsch PetscCall(DMCreateRestriction(next->dm, snes->dm, &fas->restrct)); 1240a7266b2SPatrick Farrell } else { 1259566063dSJacob Faibussowitsch PetscCall(PetscObjectReference((PetscObject)fas->interpolate)); 12679d9a41aSPeter Brune fas->restrct = fas->interpolate; 12779d9a41aSPeter Brune } 128bccf9bb3SJed Brown } 1290a7266b2SPatrick Farrell } 13079d9a41aSPeter Brune /* set the injection from the DM */ 13179d9a41aSPeter Brune if (!fas->inject) { 1329566063dSJacob Faibussowitsch PetscCall(DMHasCreateInjection(next->dm, &hasCreateInjection)); 133f833ba53SLisandro Dalcin if (hasCreateInjection) { 1349566063dSJacob Faibussowitsch PetscCall(DMCreateInjection(next->dm, snes->dm, &fas->inject)); 13579d9a41aSPeter Brune } 13679d9a41aSPeter Brune } 13779d9a41aSPeter Brune } 13823e68893SLawrence Mitchell } 139f833ba53SLisandro Dalcin 14079d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 14179d9a41aSPeter Brune if (fas->galerkin) { 1421aa26658SKarl Rupp if (next) { 1439566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next)); 1441aa26658SKarl Rupp } 1451aa26658SKarl Rupp if (fas->smoothd && fas->level != fas->levels - 1) { 1469566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes)); 1471aa26658SKarl Rupp } 1481aa26658SKarl Rupp if (fas->smoothu && fas->level != fas->levels - 1) { 1499566063dSJacob Faibussowitsch PetscCall(SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes)); 1501aa26658SKarl Rupp } 15179d9a41aSPeter Brune } 15279d9a41aSPeter Brune 153534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 154534ebe21SPeter Brune if (fas->smoothd) { 155bc3f2f05SPeter Brune if (fas->level == 0 && fas->levels != 1) { 1569566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE)); 157534ebe21SPeter Brune } else { 1589566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY)); 159534ebe21SPeter Brune } 1609566063dSJacob Faibussowitsch PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothd)); 161534ebe21SPeter Brune } 162534ebe21SPeter Brune 163534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 164534ebe21SPeter Brune if (fas->smoothu) { 165534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 1669566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE)); 167534ebe21SPeter Brune } else { 1689566063dSJacob Faibussowitsch PetscCall(SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY)); 169534ebe21SPeter Brune } 1709566063dSJacob Faibussowitsch PetscCall(SNESFASCycleSetUpSmoother_Private(snes, fas->smoothu)); 171534ebe21SPeter Brune } 172d06165b7SPeter Brune 173ab8d36c9SPeter Brune if (next) { 17479d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 1759566063dSJacob Faibussowitsch if (!next->vec_sol) PetscCall(SNESFASCreateCoarseVec(snes,&next->vec_sol)); 1769566063dSJacob Faibussowitsch if (!next->vec_rhs) PetscCall(SNESFASCreateCoarseVec(snes,&next->vec_rhs)); 1779566063dSJacob Faibussowitsch PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next)); 1789566063dSJacob Faibussowitsch PetscCall(SNESFASSetUpLineSearch_Private(snes, next)); 1799566063dSJacob Faibussowitsch PetscCall(SNESSetUp(next)); 18079d9a41aSPeter Brune } 181f833ba53SLisandro Dalcin 1826273346dSPeter Brune /* setup FAS work vectors */ 1836273346dSPeter Brune if (fas->galerkin) { 1849566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &fas->Xg)); 1859566063dSJacob Faibussowitsch PetscCall(VecDuplicate(snes->vec_sol, &fas->Fg)); 1866273346dSPeter Brune } 187421d9b32SPeter Brune PetscFunctionReturn(0); 188421d9b32SPeter Brune } 189421d9b32SPeter Brune 19040244768SBarry Smith static PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes) 191421d9b32SPeter Brune { 192ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 193ee78dd50SPeter Brune PetscInt levels = 1; 19487f44e3fSPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE; 19507144faaSPeter Brune SNESFASType fastype; 196fde0ff24SPeter Brune const char *optionsprefix; 197f1c6b773SPeter Brune SNESLineSearch linesearch; 19866585501SPeter Brune PetscInt m, n_up, n_down; 199ab8d36c9SPeter Brune SNES next; 200ab8d36c9SPeter Brune PetscBool isFine; 201421d9b32SPeter Brune 202421d9b32SPeter Brune PetscFunctionBegin; 2039566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 204d0609cedSBarry Smith PetscOptionsHeadBegin(PetscOptionsObject,"SNESFAS Options-----------------------------------"); 205ee78dd50SPeter Brune 206ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 207ab8d36c9SPeter Brune if (isFine) { 2089566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg)); 209c732cbdbSBarry Smith if (!flg && snes->dm) { 2109566063dSJacob Faibussowitsch PetscCall(DMGetRefineLevel(snes->dm,&levels)); 211c732cbdbSBarry Smith levels++; 212d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 213c732cbdbSBarry Smith } 2149566063dSJacob Faibussowitsch PetscCall(SNESFASSetLevels(snes, levels, NULL)); 21507144faaSPeter Brune fastype = fas->fastype; 2169566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg)); 21707144faaSPeter Brune if (flg) { 2189566063dSJacob Faibussowitsch PetscCall(SNESFASSetType(snes, fastype)); 21907144faaSPeter Brune } 220ee78dd50SPeter Brune 2219566063dSJacob Faibussowitsch PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix)); 2229566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg)); 223ab8d36c9SPeter Brune if (flg) { 2249566063dSJacob Faibussowitsch PetscCall(SNESFASSetCycles(snes, m)); 225fde0ff24SPeter Brune } 2269566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg)); 22787f44e3fSPeter Brune if (flg) { 2289566063dSJacob Faibussowitsch PetscCall(SNESFASSetContinuation(snes,continuationflg)); 22987f44e3fSPeter Brune } 230fde0ff24SPeter Brune 2319566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg)); 232ab8d36c9SPeter Brune if (flg) { 2339566063dSJacob Faibussowitsch PetscCall(SNESFASSetGalerkin(snes, galerkinflg)); 234ab8d36c9SPeter Brune } 235ee78dd50SPeter Brune 236928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 2379566063dSJacob 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)); 2389566063dSJacob Faibussowitsch if (flg) PetscCall(SNESFASFullSetDownSweep(snes,fas->full_downsweep)); 2399566063dSJacob 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)); 2409566063dSJacob Faibussowitsch if (flg) PetscCall(SNESFASFullSetTotal(snes,fas->full_total)); 241928e959bSPeter Brune } 242928e959bSPeter Brune 2439566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg)); 244162d76ddSPeter Brune 2459566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg)); 246162d76ddSPeter Brune 247d142ab34SLawrence Mitchell { 248d142ab34SLawrence Mitchell PetscViewer viewer; 249d142ab34SLawrence Mitchell PetscViewerFormat format; 2509566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_fas_monitor",&viewer,&format,&monflg)); 251d142ab34SLawrence Mitchell if (monflg) { 252d142ab34SLawrence Mitchell PetscViewerAndFormat *vf; 2539566063dSJacob Faibussowitsch PetscCall(PetscViewerAndFormatCreate(viewer,format,&vf)); 2549566063dSJacob Faibussowitsch PetscCall(PetscObjectDereference((PetscObject)viewer)); 2559566063dSJacob Faibussowitsch PetscCall(SNESFASSetMonitor(snes,vf,PETSC_TRUE)); 256d142ab34SLawrence Mitchell } 257d142ab34SLawrence Mitchell } 2580dd27c6cSPeter Brune flg = PETSC_FALSE; 2590dd27c6cSPeter Brune monflg = PETSC_TRUE; 2609566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg)); 2619566063dSJacob Faibussowitsch if (flg) PetscCall(SNESFASSetLog(snes,monflg)); 262ab8d36c9SPeter Brune } 263ee78dd50SPeter Brune 264d0609cedSBarry Smith PetscOptionsHeadEnd(); 265f833ba53SLisandro Dalcin 2668cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 267162d76ddSPeter Brune if (upflg) { 2689566063dSJacob Faibussowitsch PetscCall(SNESFASSetNumberSmoothUp(snes,n_up)); 269162d76ddSPeter Brune } 270162d76ddSPeter Brune if (downflg) { 2719566063dSJacob Faibussowitsch PetscCall(SNESFASSetNumberSmoothDown(snes,n_down)); 272162d76ddSPeter Brune } 273eff52c0eSPeter Brune 2749e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 2759e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 2769e764e56SPeter Brune if (!snes->linesearch) { 2779566063dSJacob Faibussowitsch PetscCall(SNESGetLineSearch(snes, &linesearch)); 2789566063dSJacob Faibussowitsch PetscCall(SNESLineSearchSetType(linesearch, SNESLINESEARCHL2)); 2799e764e56SPeter Brune } 2809e764e56SPeter Brune } 2819e764e56SPeter Brune 282ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 2839566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 2849566063dSJacob Faibussowitsch if (next) PetscCall(SNESSetFromOptions(next)); 285421d9b32SPeter Brune PetscFunctionReturn(0); 286421d9b32SPeter Brune } 287421d9b32SPeter Brune 2889804daf3SBarry Smith #include <petscdraw.h> 28940244768SBarry Smith static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 290421d9b32SPeter Brune { 291421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 292656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 293ab8d36c9SPeter Brune PetscInt i; 294ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 295421d9b32SPeter Brune 296421d9b32SPeter Brune PetscFunctionBegin; 2979566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 298ab8d36c9SPeter Brune if (isFine) { 2999566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii)); 3009566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw)); 301421d9b32SPeter Brune if (iascii) { 302*63a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, " type is %s, levels=%" PetscInt_FMT ", cycles=%" PetscInt_FMT "\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles)); 303ab8d36c9SPeter Brune if (fas->galerkin) { 3049566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n")); 305421d9b32SPeter Brune } else { 3069566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n")); 307421d9b32SPeter Brune } 308ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 3099566063dSJacob Faibussowitsch PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes)); 3109566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetSmootherUp(levelsnes, &smoothu)); 3119566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetSmootherDown(levelsnes, &smoothd)); 312ab8d36c9SPeter Brune if (!i) { 313*63a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Coarse grid solver -- level %" PetscInt_FMT " -------------------------------\n",i)); 314421d9b32SPeter Brune } else { 315*63a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Down solver (pre-smoother) on level %" PetscInt_FMT " -------------------------------\n",i)); 316421d9b32SPeter Brune } 3179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 318166b3ea4SJed Brown if (smoothd) { 3199566063dSJacob Faibussowitsch PetscCall(SNESView(smoothd,viewer)); 320166b3ea4SJed Brown } else { 3219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Not yet available\n")); 322166b3ea4SJed Brown } 3239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 324ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 3259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Up solver (post-smoother) same as down solver (pre-smoother)\n")); 326ab8d36c9SPeter Brune } else if (i) { 327*63a3b9bcSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer," Up solver (post-smoother) on level %" PetscInt_FMT " -------------------------------\n",i)); 3289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 329166b3ea4SJed Brown if (smoothu) { 3309566063dSJacob Faibussowitsch PetscCall(SNESView(smoothu,viewer)); 331166b3ea4SJed Brown } else { 3329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer,"Not yet available\n")); 333166b3ea4SJed Brown } 3349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 335ab8d36c9SPeter Brune } 336ab8d36c9SPeter Brune } 337656ede7eSPeter Brune } else if (isdraw) { 338656ede7eSPeter Brune PetscDraw draw; 339b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 340656ede7eSPeter Brune SNES_FAS *curfas = fas; 3419566063dSJacob Faibussowitsch PetscCall(PetscViewerDrawGetDraw(viewer,0,&draw)); 3429566063dSJacob Faibussowitsch PetscCall(PetscDrawGetCurrentPoint(draw,&x,&y)); 3439566063dSJacob Faibussowitsch PetscCall(PetscDrawStringGetSize(draw,&wth,&th)); 344656ede7eSPeter Brune bottom = y - th; 345656ede7eSPeter Brune while (curfas) { 346b4375e8dSPeter Brune if (!curfas->smoothu) { 3479566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw,x,bottom)); 3489566063dSJacob Faibussowitsch if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd,viewer)); 3499566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 350b4375e8dSPeter Brune } else { 351b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 3529566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw,x-w,bottom)); 3539566063dSJacob Faibussowitsch if (curfas->smoothd) PetscCall(SNESView(curfas->smoothd,viewer)); 3549566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 3559566063dSJacob Faibussowitsch PetscCall(PetscDrawPushCurrentPoint(draw,x+w,bottom)); 3569566063dSJacob Faibussowitsch if (curfas->smoothu) PetscCall(SNESView(curfas->smoothu,viewer)); 3579566063dSJacob Faibussowitsch PetscCall(PetscDrawPopCurrentPoint(draw)); 358b4375e8dSPeter Brune } 359656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 360656ede7eSPeter Brune bottom -= 5*th; 3611aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 3620298fd71SBarry Smith else curfas = NULL; 363656ede7eSPeter Brune } 364421d9b32SPeter Brune } 365ab8d36c9SPeter Brune } 366421d9b32SPeter Brune PetscFunctionReturn(0); 367421d9b32SPeter Brune } 368421d9b32SPeter Brune 36939bd7f45SPeter Brune /* 37039bd7f45SPeter Brune Defines the action of the downsmoother 37139bd7f45SPeter Brune */ 37240244768SBarry Smith static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 373b9c2fdf1SPeter Brune { 374742fe5e2SPeter Brune SNESConvergedReason reason; 375ab8d36c9SPeter Brune Vec FPC; 376ab8d36c9SPeter Brune SNES smoothd; 3776cbb2f26SLawrence Mitchell PetscBool flg; 3780dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 3796e111a19SKarl Rupp 380421d9b32SPeter Brune PetscFunctionBegin; 3819566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetSmootherDown(snes, &smoothd)); 3829566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(smoothd, F)); 3839566063dSJacob Faibussowitsch if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve,smoothd,B,X,0)); 3849566063dSJacob Faibussowitsch PetscCall(SNESSolve(smoothd, B, X)); 3859566063dSJacob Faibussowitsch if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve,smoothd,B,X,0)); 386742fe5e2SPeter Brune /* check convergence reason for the smoother */ 3879566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(smoothd,&reason)); 3881ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 389742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 390742fe5e2SPeter Brune PetscFunctionReturn(0); 391742fe5e2SPeter Brune } 3926cbb2f26SLawrence Mitchell 3939566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(smoothd, &FPC, NULL, NULL)); 3949566063dSJacob Faibussowitsch PetscCall(SNESGetAlwaysComputesFinalResidual(smoothd, &flg)); 3956cbb2f26SLawrence Mitchell if (!flg) { 3969566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(smoothd, X, FPC)); 3976cbb2f26SLawrence Mitchell } 3989566063dSJacob Faibussowitsch PetscCall(VecCopy(FPC, F)); 3999566063dSJacob Faibussowitsch if (fnorm) PetscCall(VecNorm(F,NORM_2,fnorm)); 40039bd7f45SPeter Brune PetscFunctionReturn(0); 40139bd7f45SPeter Brune } 40239bd7f45SPeter Brune 40339bd7f45SPeter Brune /* 40407144faaSPeter Brune Defines the action of the upsmoother 40539bd7f45SPeter Brune */ 40640244768SBarry Smith static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 4070adebc6cSBarry Smith { 40839bd7f45SPeter Brune SNESConvergedReason reason; 409ab8d36c9SPeter Brune Vec FPC; 410ab8d36c9SPeter Brune SNES smoothu; 4116cbb2f26SLawrence Mitchell PetscBool flg; 4120dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 413ab8d36c9SPeter Brune 4146e111a19SKarl Rupp PetscFunctionBegin; 4159566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetSmootherUp(snes, &smoothu)); 4169566063dSJacob Faibussowitsch if (fas->eventsmoothsolve) PetscCall(PetscLogEventBegin(fas->eventsmoothsolve,smoothu,0,0,0)); 4179566063dSJacob Faibussowitsch PetscCall(SNESSolve(smoothu, B, X)); 4189566063dSJacob Faibussowitsch if (fas->eventsmoothsolve) PetscCall(PetscLogEventEnd(fas->eventsmoothsolve,smoothu,0,0,0)); 41939bd7f45SPeter Brune /* check convergence reason for the smoother */ 4209566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(smoothu,&reason)); 4211ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 42239bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 42339bd7f45SPeter Brune PetscFunctionReturn(0); 42439bd7f45SPeter Brune } 4259566063dSJacob Faibussowitsch PetscCall(SNESGetFunction(smoothu, &FPC, NULL, NULL)); 4269566063dSJacob Faibussowitsch PetscCall(SNESGetAlwaysComputesFinalResidual(smoothu, &flg)); 4276cbb2f26SLawrence Mitchell if (!flg) { 4289566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(smoothu, X, FPC)); 4296cbb2f26SLawrence Mitchell } 4309566063dSJacob Faibussowitsch PetscCall(VecCopy(FPC, F)); 4319566063dSJacob Faibussowitsch if (fnorm) PetscCall(VecNorm(F,NORM_2,fnorm)); 43239bd7f45SPeter Brune PetscFunctionReturn(0); 43339bd7f45SPeter Brune } 43439bd7f45SPeter Brune 435938e4a01SJed Brown /*@ 436938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 437938e4a01SJed Brown 438938e4a01SJed Brown Collective 439938e4a01SJed Brown 4404165533cSJose E. Roman Input Parameter: 441938e4a01SJed Brown . snes - SNESFAS 442938e4a01SJed Brown 4434165533cSJose E. Roman Output Parameter: 444938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 445938e4a01SJed Brown 446938e4a01SJed Brown Level: developer 447938e4a01SJed Brown 448938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 449938e4a01SJed Brown @*/ 450938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 451938e4a01SJed Brown { 452f833ba53SLisandro Dalcin SNES_FAS *fas; 453938e4a01SJed Brown 454938e4a01SJed Brown PetscFunctionBegin; 455f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 456f833ba53SLisandro Dalcin PetscValidPointer(Xcoarse,2); 457f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 4581aa26658SKarl Rupp if (fas->rscale) { 4599566063dSJacob Faibussowitsch PetscCall(VecDuplicate(fas->rscale,Xcoarse)); 460f5af7f23SKarl Rupp } else if (fas->inject) { 4619566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(fas->inject,Xcoarse,NULL)); 46213903a91SSatish Balay } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection"); 463938e4a01SJed Brown PetscFunctionReturn(0); 464938e4a01SJed Brown } 465938e4a01SJed Brown 466e9923e8dSJed Brown /*@ 467e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 468e9923e8dSJed Brown 469e9923e8dSJed Brown Collective 470e9923e8dSJed Brown 4714165533cSJose E. Roman Input Parameters: 472e9923e8dSJed Brown + fine - SNES from which to restrict 473e9923e8dSJed Brown - Xfine - vector to restrict 474e9923e8dSJed Brown 4754165533cSJose E. Roman Output Parameter: 476e9923e8dSJed Brown . Xcoarse - result of restriction 477e9923e8dSJed Brown 478e9923e8dSJed Brown Level: developer 479e9923e8dSJed Brown 480e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 481e9923e8dSJed Brown @*/ 482e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 483e9923e8dSJed Brown { 484f833ba53SLisandro Dalcin SNES_FAS *fas; 485e9923e8dSJed Brown 486e9923e8dSJed Brown PetscFunctionBegin; 487f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(fine,SNES_CLASSID,1,SNESFAS); 488e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 489e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 490f833ba53SLisandro Dalcin fas = (SNES_FAS*)fine->data; 491e9923e8dSJed Brown if (fas->inject) { 4929566063dSJacob Faibussowitsch PetscCall(MatRestrict(fas->inject,Xfine,Xcoarse)); 493e9923e8dSJed Brown } else { 4949566063dSJacob Faibussowitsch PetscCall(MatRestrict(fas->restrct,Xfine,Xcoarse)); 4959566063dSJacob Faibussowitsch PetscCall(VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse)); 496e9923e8dSJed Brown } 497e9923e8dSJed Brown PetscFunctionReturn(0); 498e9923e8dSJed Brown } 499e9923e8dSJed Brown 50039bd7f45SPeter Brune /* 50139bd7f45SPeter Brune 5026dfbafefSToby Isaac Performs a variant of FAS using the interpolated total coarse solution 5036dfbafefSToby Isaac 5046dfbafefSToby Isaac fine problem: F(x) = b 5056dfbafefSToby Isaac coarse problem: F^c(x^c) = Rb, Initial guess Rx 5066dfbafefSToby Isaac interpolated solution: x^f = I x^c (total solution interpolation 5076dfbafefSToby Isaac 5086dfbafefSToby Isaac */ 5096dfbafefSToby Isaac static PetscErrorCode SNESFASInterpolatedCoarseSolution(SNES snes, Vec X, Vec X_new) 5106dfbafefSToby Isaac { 5116dfbafefSToby Isaac Vec X_c, B_c; 5126dfbafefSToby Isaac SNESConvergedReason reason; 5136dfbafefSToby Isaac SNES next; 5146dfbafefSToby Isaac Mat restrct, interpolate; 5156dfbafefSToby Isaac SNES_FAS *fasc; 5166dfbafefSToby Isaac 5176dfbafefSToby Isaac PetscFunctionBegin; 5189566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 5196dfbafefSToby Isaac if (next) { 5206dfbafefSToby Isaac fasc = (SNES_FAS*)next->data; 5216dfbafefSToby Isaac 5229566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetRestriction(snes, &restrct)); 5239566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate)); 5246dfbafefSToby Isaac 5256dfbafefSToby Isaac X_c = next->vec_sol; 5266dfbafefSToby Isaac 5279566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0)); 5286dfbafefSToby Isaac /* restrict the total solution: Rb */ 5299566063dSJacob Faibussowitsch PetscCall(SNESFASRestrict(snes, X, X_c)); 5306dfbafefSToby Isaac B_c = next->vec_rhs; 5316dfbafefSToby Isaac if (snes->vec_rhs) { 5326dfbafefSToby Isaac /* restrict the total rhs defect: Rb */ 5339566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, snes->vec_rhs, B_c)); 5346dfbafefSToby Isaac } else { 5359566063dSJacob Faibussowitsch PetscCall(VecSet(B_c, 0.)); 5366dfbafefSToby Isaac } 5379566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0)); 5386dfbafefSToby Isaac 5399566063dSJacob Faibussowitsch PetscCall(SNESSolve(next, B_c, X_c)); 5409566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next,&reason)); 5416dfbafefSToby Isaac if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 5426dfbafefSToby Isaac snes->reason = SNES_DIVERGED_INNER; 5436dfbafefSToby Isaac PetscFunctionReturn(0); 5446dfbafefSToby Isaac } 5456dfbafefSToby Isaac /* x^f <- Ix^c*/ 5466dfbafefSToby Isaac DM dmc,dmf; 5476dfbafefSToby Isaac 5489566063dSJacob Faibussowitsch PetscCall(SNESGetDM(next, &dmc)); 5499566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &dmf)); 5509566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0)); 5519566063dSJacob Faibussowitsch PetscCall(DMInterpolateSolution(dmc, dmf, interpolate, X_c, X_new)); 5529566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0)); 5539566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) X_c, "Coarse solution")); 5549566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view")); 5559566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) X_new, "Updated Fine solution")); 5569566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view")); 5576dfbafefSToby Isaac } 5586dfbafefSToby Isaac PetscFunctionReturn(0); 5596dfbafefSToby Isaac } 5606dfbafefSToby Isaac 5616dfbafefSToby Isaac /* 5626dfbafefSToby Isaac 56339bd7f45SPeter Brune Performs the FAS coarse correction as: 56439bd7f45SPeter Brune 565b5527d98SMatthew G. Knepley fine problem: F(x) = b 566b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c 56739bd7f45SPeter Brune 568b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 56939bd7f45SPeter Brune 57039bd7f45SPeter Brune */ 5710adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 5720adebc6cSBarry Smith { 57339bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 57439bd7f45SPeter Brune SNESConvergedReason reason; 575ab8d36c9SPeter Brune SNES next; 576ab8d36c9SPeter Brune Mat restrct, interpolate; 5770dd27c6cSPeter Brune SNES_FAS *fasc; 5785fd66863SKarl Rupp 57939bd7f45SPeter Brune PetscFunctionBegin; 5809566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 581ab8d36c9SPeter Brune if (next) { 5820dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 5830dd27c6cSPeter Brune 5849566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetRestriction(snes, &restrct)); 5859566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate)); 586ab8d36c9SPeter Brune 587ab8d36c9SPeter Brune X_c = next->vec_sol; 588ab8d36c9SPeter Brune Xo_c = next->work[0]; 589ab8d36c9SPeter Brune F_c = next->vec_func; 590ab8d36c9SPeter Brune B_c = next->vec_rhs; 591efe1f98aSPeter Brune 5929566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0)); 5939566063dSJacob Faibussowitsch PetscCall(SNESFASRestrict(snes, X, Xo_c)); 5945609cbf2SMatthew G. Knepley /* restrict the defect: R(F(x) - b) */ 5959566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, F, B_c)); 5969566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0)); 5970dd27c6cSPeter Brune 5989566063dSJacob Faibussowitsch if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual,next,0,0,0)); 5995609cbf2SMatthew G. Knepley /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */ 6009566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(next, Xo_c, F_c)); 6019566063dSJacob Faibussowitsch if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual,next,0,0,0)); 6020dd27c6cSPeter Brune 6030dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 6049566063dSJacob Faibussowitsch PetscCall(VecCopy(B_c, X_c)); 6059566063dSJacob Faibussowitsch PetscCall(VecCopy(F_c, B_c)); 6069566063dSJacob Faibussowitsch PetscCall(VecCopy(X_c, F_c)); 607ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 6089566063dSJacob Faibussowitsch PetscCall(VecCopy(Xo_c, X_c)); 609c90fad12SPeter Brune 610c90fad12SPeter Brune /* recurse to the next level */ 6119566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next, F_c)); 6129566063dSJacob Faibussowitsch PetscCall(SNESSolve(next, B_c, X_c)); 6139566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next,&reason)); 614742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 615742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 616742fe5e2SPeter Brune PetscFunctionReturn(0); 617742fe5e2SPeter Brune } 618fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 6199566063dSJacob Faibussowitsch PetscCall(VecAXPY(X_c, -1.0, Xo_c)); 6200dd27c6cSPeter Brune 6219566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0)); 6229566063dSJacob Faibussowitsch PetscCall(MatInterpolateAdd(interpolate, X_c, X, X_new)); 6239566063dSJacob Faibussowitsch if (fasc->eventinterprestrict) PetscCall(PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0)); 6249566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) X_c, "Coarse correction")); 6259566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view")); 6269566063dSJacob Faibussowitsch PetscCall(PetscObjectSetName((PetscObject) X_new, "Updated Fine solution")); 6279566063dSJacob Faibussowitsch PetscCall(VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view")); 628293a7e31SPeter Brune } 62939bd7f45SPeter Brune PetscFunctionReturn(0); 63039bd7f45SPeter Brune } 63139bd7f45SPeter Brune 63239bd7f45SPeter Brune /* 63339bd7f45SPeter Brune 63439bd7f45SPeter Brune The additive cycle looks like: 63539bd7f45SPeter Brune 63607144faaSPeter Brune xhat = x 63707144faaSPeter Brune xhat = dS(x, b) 63807144faaSPeter Brune x = coarsecorrection(xhat, b_d) 63907144faaSPeter Brune x = x + nu*(xhat - x); 64039bd7f45SPeter Brune (optional) x = uS(x, b) 64139bd7f45SPeter Brune 64239bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 64339bd7f45SPeter Brune 64439bd7f45SPeter Brune */ 64540244768SBarry Smith static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 6460adebc6cSBarry Smith { 64707144faaSPeter Brune Vec F, B, Xhat; 64822c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 64907144faaSPeter Brune SNESConvergedReason reason; 65022c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 651422a814eSBarry Smith SNESLineSearchReason lsresult; 652ab8d36c9SPeter Brune SNES next; 653ab8d36c9SPeter Brune Mat restrct, interpolate; 6540dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 6550dd27c6cSPeter Brune 65639bd7f45SPeter Brune PetscFunctionBegin; 6579566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 65839bd7f45SPeter Brune F = snes->vec_func; 65939bd7f45SPeter Brune B = snes->vec_rhs; 660e7f468e7SPeter Brune Xhat = snes->work[1]; 6619566063dSJacob Faibussowitsch PetscCall(VecCopy(X, Xhat)); 66207144faaSPeter Brune /* recurse first */ 663ab8d36c9SPeter Brune if (next) { 6640dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 6659566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetRestriction(snes, &restrct)); 6669566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetInterpolation(snes, &interpolate)); 6679566063dSJacob Faibussowitsch if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual,snes,0,0,0)); 6689566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes, Xhat, F)); 6699566063dSJacob Faibussowitsch if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual,snes,0,0,0)); 6709566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); 671ab8d36c9SPeter Brune X_c = next->vec_sol; 672ab8d36c9SPeter Brune Xo_c = next->work[0]; 673ab8d36c9SPeter Brune F_c = next->vec_func; 674ab8d36c9SPeter Brune B_c = next->vec_rhs; 67539bd7f45SPeter Brune 6769566063dSJacob Faibussowitsch PetscCall(SNESFASRestrict(snes,Xhat,Xo_c)); 67707144faaSPeter Brune /* restrict the defect */ 6789566063dSJacob Faibussowitsch PetscCall(MatRestrict(restrct, F, B_c)); 67907144faaSPeter Brune 68007144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 6819566063dSJacob Faibussowitsch if (fasc->eventresidual) PetscCall(PetscLogEventBegin(fasc->eventresidual,next,0,0,0)); 6829566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(next, Xo_c, F_c)); 6839566063dSJacob Faibussowitsch if (fasc->eventresidual) PetscCall(PetscLogEventEnd(fasc->eventresidual,next,0,0,0)); 6849566063dSJacob Faibussowitsch PetscCall(VecCopy(B_c, X_c)); 6859566063dSJacob Faibussowitsch PetscCall(VecCopy(F_c, B_c)); 6869566063dSJacob Faibussowitsch PetscCall(VecCopy(X_c, F_c)); 68707144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 6889566063dSJacob Faibussowitsch PetscCall(VecCopy(Xo_c, X_c)); 68907144faaSPeter Brune 69007144faaSPeter Brune /* recurse */ 6919566063dSJacob Faibussowitsch PetscCall(SNESSetInitialFunction(next, F_c)); 6929566063dSJacob Faibussowitsch PetscCall(SNESSolve(next, B_c, X_c)); 69307144faaSPeter Brune 69407144faaSPeter Brune /* smooth on this level */ 6959566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &fnorm)); 69607144faaSPeter Brune 6979566063dSJacob Faibussowitsch PetscCall(SNESGetConvergedReason(next,&reason)); 69807144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 69907144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 70007144faaSPeter Brune PetscFunctionReturn(0); 70107144faaSPeter Brune } 70207144faaSPeter Brune 70307144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 7049566063dSJacob Faibussowitsch PetscCall(VecAYPX(X_c, -1.0, Xo_c)); 7059566063dSJacob Faibussowitsch PetscCall(MatInterpolate(interpolate, X_c, Xhat)); 70607144faaSPeter Brune 707ddebd997SPeter Brune /* additive correction of the coarse direction*/ 7089566063dSJacob Faibussowitsch PetscCall(SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat)); 7099566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetReason(snes->linesearch, &lsresult)); 7109566063dSJacob Faibussowitsch PetscCall(SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm)); 711422a814eSBarry Smith if (lsresult) { 7129e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 7139e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 7149e764e56SPeter Brune PetscFunctionReturn(0); 7159e764e56SPeter Brune } 7169e764e56SPeter Brune } 71707144faaSPeter Brune } else { 7189566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 71907144faaSPeter Brune } 72039bd7f45SPeter Brune PetscFunctionReturn(0); 72139bd7f45SPeter Brune } 72239bd7f45SPeter Brune 72339bd7f45SPeter Brune /* 72439bd7f45SPeter Brune 72539bd7f45SPeter Brune Defines the FAS cycle as: 72639bd7f45SPeter Brune 727b5527d98SMatthew G. Knepley fine problem: F(x) = b 72839bd7f45SPeter Brune coarse problem: F^c(x) = b^c 72939bd7f45SPeter Brune 730b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 73139bd7f45SPeter Brune 73239bd7f45SPeter Brune correction: 73339bd7f45SPeter Brune 73439bd7f45SPeter Brune x = x + I(x^c - Rx) 73539bd7f45SPeter Brune 73639bd7f45SPeter Brune */ 73740244768SBarry Smith static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7380adebc6cSBarry Smith { 73939bd7f45SPeter Brune 74039bd7f45SPeter Brune Vec F,B; 74134d65b3cSPeter Brune SNES next; 74239bd7f45SPeter Brune 74339bd7f45SPeter Brune PetscFunctionBegin; 74439bd7f45SPeter Brune F = snes->vec_func; 74539bd7f45SPeter Brune B = snes->vec_rhs; 74639bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 7479566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes, &next)); 7489566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm)); 74934d65b3cSPeter Brune if (next) { 7509566063dSJacob Faibussowitsch PetscCall(SNESFASCoarseCorrection(snes, X, F, X)); 7519566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm)); 752fe6f9142SPeter Brune } 753fa9694d7SPeter Brune PetscFunctionReturn(0); 754421d9b32SPeter Brune } 755421d9b32SPeter Brune 75640244768SBarry Smith static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 7578c94862eSPeter Brune { 7588c94862eSPeter Brune SNES next; 7598c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 7608c94862eSPeter Brune PetscBool isFine; 7618c94862eSPeter Brune 7628c94862eSPeter Brune PetscFunctionBegin; 7638c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 7649566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes,&isFine)); 7659566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes,&next)); 7668c94862eSPeter Brune fas->full_stage = 0; 7679566063dSJacob Faibussowitsch if (next) PetscCall(SNESFASCycleSetupPhase_Full(next)); 7688c94862eSPeter Brune PetscFunctionReturn(0); 7698c94862eSPeter Brune } 7708c94862eSPeter Brune 77140244768SBarry Smith static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 7728c94862eSPeter Brune { 7738c94862eSPeter Brune Vec F,B; 7748c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 7758c94862eSPeter Brune PetscBool isFine; 7768c94862eSPeter Brune SNES next; 7778c94862eSPeter Brune 7788c94862eSPeter Brune PetscFunctionBegin; 7798c94862eSPeter Brune F = snes->vec_func; 7808c94862eSPeter Brune B = snes->vec_rhs; 7819566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes,&isFine)); 7829566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes,&next)); 7838c94862eSPeter Brune 7848c94862eSPeter Brune if (isFine) { 7859566063dSJacob Faibussowitsch PetscCall(SNESFASCycleSetupPhase_Full(snes)); 7868c94862eSPeter Brune } 7878c94862eSPeter Brune 7888c94862eSPeter Brune if (fas->full_stage == 0) { 789928e959bSPeter Brune /* downsweep */ 7908c94862eSPeter Brune if (next) { 7918c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 7929566063dSJacob Faibussowitsch if (fas->full_downsweep) PetscCall(SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm)); 793e140b65aSPatrick Farrell fas->full_downsweep = PETSC_TRUE; 7949566063dSJacob Faibussowitsch if (fas->full_total) PetscCall(SNESFASInterpolatedCoarseSolution(snes,X,X)); 7959566063dSJacob Faibussowitsch else PetscCall(SNESFASCoarseCorrection(snes,X,F,X)); 7966dfbafefSToby Isaac fas->full_total = PETSC_FALSE; 7979566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm)); 7988c94862eSPeter Brune if (fas->level != 1) next->max_its -= 1; 7998c94862eSPeter Brune } else { 800a3a80b83SMatthew G. Knepley /* The smoother on the coarse level is the coarse solver */ 8019566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm)); 8028c94862eSPeter Brune } 8038c94862eSPeter Brune fas->full_stage = 1; 8048c94862eSPeter Brune } else if (fas->full_stage == 1) { 8059566063dSJacob Faibussowitsch if (snes->iter == 0) PetscCall(SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm)); 8068c94862eSPeter Brune if (next) { 8079566063dSJacob Faibussowitsch PetscCall(SNESFASCoarseCorrection(snes,X,F,X)); 8089566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm)); 8098c94862eSPeter Brune } 8108c94862eSPeter Brune } 8118c94862eSPeter Brune /* final v-cycle */ 8128c94862eSPeter Brune if (isFine) { 8138c94862eSPeter Brune if (next) { 8149566063dSJacob Faibussowitsch PetscCall(SNESFASCoarseCorrection(snes,X,F,X)); 8159566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm)); 8168c94862eSPeter Brune } 8178c94862eSPeter Brune } 8188c94862eSPeter Brune PetscFunctionReturn(0); 8198c94862eSPeter Brune } 8208c94862eSPeter Brune 82140244768SBarry Smith static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 82234d65b3cSPeter Brune { 82334d65b3cSPeter Brune Vec F,B; 82434d65b3cSPeter Brune SNES next; 82534d65b3cSPeter Brune 82634d65b3cSPeter Brune PetscFunctionBegin; 82734d65b3cSPeter Brune F = snes->vec_func; 82834d65b3cSPeter Brune B = snes->vec_rhs; 8299566063dSJacob Faibussowitsch PetscCall(SNESFASCycleGetCorrection(snes,&next)); 83034d65b3cSPeter Brune if (next) { 8319566063dSJacob Faibussowitsch PetscCall(SNESFASCoarseCorrection(snes,X,F,X)); 8329566063dSJacob Faibussowitsch PetscCall(SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm)); 83334d65b3cSPeter Brune } else { 8349566063dSJacob Faibussowitsch PetscCall(SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm)); 83534d65b3cSPeter Brune } 83634d65b3cSPeter Brune PetscFunctionReturn(0); 83734d65b3cSPeter Brune } 83834d65b3cSPeter Brune 839fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE; 840fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n" 841fffbeea8SBarry Smith " title = {Composing Scalable Nonlinear Algebraic Solvers},\n" 842fffbeea8SBarry Smith " author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n" 843fffbeea8SBarry Smith " year = 2013,\n" 844fffbeea8SBarry Smith " type = Preprint,\n" 845fffbeea8SBarry Smith " number = {ANL/MCS-P2010-0112},\n" 846fffbeea8SBarry Smith " institution = {Argonne National Laboratory}\n}\n"; 847fffbeea8SBarry Smith 84840244768SBarry Smith static PetscErrorCode SNESSolve_FAS(SNES snes) 849421d9b32SPeter Brune { 850f833ba53SLisandro Dalcin PetscInt i; 851ddb5aff1SPeter Brune Vec X, F; 852fe6f9142SPeter Brune PetscReal fnorm; 853b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 854b17ce1afSJed Brown DM dm; 855e70c42e5SPeter Brune PetscBool isFine; 856b17ce1afSJed Brown 857421d9b32SPeter Brune PetscFunctionBegin; 8582c71b3e2SJacob Faibussowitsch PetscCheckFalse(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); 859c579b300SPatrick Farrell 8609566063dSJacob Faibussowitsch PetscCall(PetscCitationsRegister(SNESCitation,&SNEScite)); 861fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 862fa9694d7SPeter Brune X = snes->vec_sol; 863f5a6d4f9SBarry Smith F = snes->vec_func; 864293a7e31SPeter Brune 8659566063dSJacob Faibussowitsch PetscCall(SNESFASCycleIsFine(snes, &isFine)); 866293a7e31SPeter Brune /* norm setup */ 8679566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 868fe6f9142SPeter Brune snes->iter = 0; 869f833ba53SLisandro Dalcin snes->norm = 0; 8709566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 871e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 8729566063dSJacob Faibussowitsch if (fas->eventresidual) PetscCall(PetscLogEventBegin(fas->eventresidual,snes,0,0,0)); 8739566063dSJacob Faibussowitsch PetscCall(SNESComputeFunction(snes,X,F)); 8749566063dSJacob Faibussowitsch if (fas->eventresidual) PetscCall(PetscLogEventEnd(fas->eventresidual,snes,0,0,0)); 8751aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 876e4ed7901SPeter Brune 8779566063dSJacob Faibussowitsch PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F|| */ 878422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 8799566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 880fe6f9142SPeter Brune snes->norm = fnorm; 8819566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 8829566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,fnorm,0)); 8839566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,snes->iter,fnorm)); 884fe6f9142SPeter Brune 885fe6f9142SPeter Brune /* test convergence */ 8869566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP)); 887fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 888e4ed7901SPeter Brune 889b9c2fdf1SPeter Brune if (isFine) { 890b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 8919566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes,&dm)); 892b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 893b17ce1afSJed Brown DM dmcoarse; 8949566063dSJacob Faibussowitsch PetscCall(SNESGetDM(ffas->next,&dmcoarse)); 8959566063dSJacob Faibussowitsch PetscCall(DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse)); 896b17ce1afSJed Brown dm = dmcoarse; 897b17ce1afSJed Brown } 898b9c2fdf1SPeter Brune } 899b17ce1afSJed Brown 900f833ba53SLisandro Dalcin for (i = 0; i < snes->max_its; i++) { 901fe6f9142SPeter Brune /* Call general purpose update function */ 902fe6f9142SPeter Brune if (snes->ops->update) { 9039566063dSJacob Faibussowitsch PetscCall((*snes->ops->update)(snes, snes->iter)); 904fe6f9142SPeter Brune } 905f833ba53SLisandro Dalcin 90607144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 9079566063dSJacob Faibussowitsch PetscCall(SNESFASCycle_Multiplicative(snes, X)); 9088c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_ADDITIVE) { 9099566063dSJacob Faibussowitsch PetscCall(SNESFASCycle_Additive(snes, X)); 9108c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_FULL) { 9119566063dSJacob Faibussowitsch PetscCall(SNESFASCycle_Full(snes, X)); 91234d65b3cSPeter Brune } else if (fas->fastype == SNES_FAS_KASKADE) { 9139566063dSJacob Faibussowitsch PetscCall(SNESFASCycle_Kaskade(snes, X)); 9146c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type"); 915742fe5e2SPeter Brune 916742fe5e2SPeter Brune /* check for FAS cycle divergence */ 9171aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 918b9c2fdf1SPeter Brune 919c90fad12SPeter Brune /* Monitor convergence */ 9209566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes)); 921c90fad12SPeter Brune snes->iter = i+1; 9229566063dSJacob Faibussowitsch PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes)); 9239566063dSJacob Faibussowitsch PetscCall(SNESLogConvergenceHistory(snes,snes->norm,0)); 9249566063dSJacob Faibussowitsch PetscCall(SNESMonitor(snes,snes->iter,snes->norm)); 925c90fad12SPeter Brune /* Test for convergence */ 92666585501SPeter Brune if (isFine) { 9279566063dSJacob Faibussowitsch PetscCall((*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP)); 928c90fad12SPeter Brune if (snes->reason) break; 929fe6f9142SPeter Brune } 93066585501SPeter Brune } 931f833ba53SLisandro Dalcin if (i == snes->max_its) { 932*63a3b9bcSJacob Faibussowitsch PetscCall(PetscInfo(snes, "Maximum number of iterations has been reached: %" PetscInt_FMT "\n", i)); 933fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 934fe6f9142SPeter Brune } 935421d9b32SPeter Brune PetscFunctionReturn(0); 936421d9b32SPeter Brune } 93740244768SBarry Smith 93840244768SBarry Smith /*MC 93940244768SBarry Smith 94040244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 94140244768SBarry Smith 94240244768SBarry Smith The nonlinear problem is solved by correction using coarse versions 94340244768SBarry Smith of the nonlinear problem. This problem is perturbed so that a projected 94440244768SBarry Smith solution of the fine problem elicits no correction from the coarse problem. 94540244768SBarry Smith 94640244768SBarry Smith Options Database: 94740244768SBarry Smith + -snes_fas_levels - The number of levels 94840244768SBarry Smith . -snes_fas_cycles<1> - The number of cycles -- 1 for V, 2 for W 94940244768SBarry Smith . -snes_fas_type<additive,multiplicative,full,kaskade> - Additive or multiplicative cycle 95040244768SBarry Smith . -snes_fas_galerkin<PETSC_FALSE> - Form coarse problems by projection back upon the fine problem 95140244768SBarry Smith . -snes_fas_smoothup<1> - The number of iterations of the post-smoother 95240244768SBarry Smith . -snes_fas_smoothdown<1> - The number of iterations of the pre-smoother 95340244768SBarry Smith . -snes_fas_monitor - Monitor progress of all of the levels 95440244768SBarry Smith . -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS 95540244768SBarry Smith . -fas_levels_snes_ - SNES options for all smoothers 95640244768SBarry Smith . -fas_levels_cycle_snes_ - SNES options for all cycles 95740244768SBarry Smith . -fas_levels_i_snes_ - SNES options for the smoothers on level i 95840244768SBarry Smith . -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i 95940244768SBarry Smith - -fas_coarse_snes_ - SNES options for the coarsest smoother 96040244768SBarry Smith 96140244768SBarry Smith Notes: 96240244768SBarry Smith The organization of the FAS solver is slightly different from the organization of PCMG 96340244768SBarry Smith As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance. 96440244768SBarry Smith The cycle SNES instance may be used for monitoring convergence on a particular level. 96540244768SBarry Smith 96640244768SBarry Smith Level: beginner 96740244768SBarry Smith 96840244768SBarry Smith References: 969606c0280SSatish Balay . * - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 97040244768SBarry Smith SIAM Review, 57(4), 2015 97140244768SBarry Smith 97240244768SBarry Smith .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 97340244768SBarry Smith M*/ 97440244768SBarry Smith 97540244768SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes) 97640244768SBarry Smith { 97740244768SBarry Smith SNES_FAS *fas; 97840244768SBarry Smith 97940244768SBarry Smith PetscFunctionBegin; 98040244768SBarry Smith snes->ops->destroy = SNESDestroy_FAS; 98140244768SBarry Smith snes->ops->setup = SNESSetUp_FAS; 98240244768SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_FAS; 98340244768SBarry Smith snes->ops->view = SNESView_FAS; 98440244768SBarry Smith snes->ops->solve = SNESSolve_FAS; 98540244768SBarry Smith snes->ops->reset = SNESReset_FAS; 98640244768SBarry Smith 98740244768SBarry Smith snes->usesksp = PETSC_FALSE; 988efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 98940244768SBarry Smith 99040244768SBarry Smith if (!snes->tolerancesset) { 99140244768SBarry Smith snes->max_funcs = 30000; 99240244768SBarry Smith snes->max_its = 10000; 99340244768SBarry Smith } 99440244768SBarry Smith 9954fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 9964fc747eaSLawrence Mitchell 9979566063dSJacob Faibussowitsch PetscCall(PetscNewLog(snes,&fas)); 99840244768SBarry Smith 99940244768SBarry Smith snes->data = (void*) fas; 100040244768SBarry Smith fas->level = 0; 100140244768SBarry Smith fas->levels = 1; 100240244768SBarry Smith fas->n_cycles = 1; 100340244768SBarry Smith fas->max_up_it = 1; 100440244768SBarry Smith fas->max_down_it = 1; 100540244768SBarry Smith fas->smoothu = NULL; 100640244768SBarry Smith fas->smoothd = NULL; 100740244768SBarry Smith fas->next = NULL; 100840244768SBarry Smith fas->previous = NULL; 100940244768SBarry Smith fas->fine = snes; 101040244768SBarry Smith fas->interpolate = NULL; 101140244768SBarry Smith fas->restrct = NULL; 101240244768SBarry Smith fas->inject = NULL; 101340244768SBarry Smith fas->usedmfornumberoflevels = PETSC_FALSE; 101440244768SBarry Smith fas->fastype = SNES_FAS_MULTIPLICATIVE; 101540244768SBarry Smith fas->full_downsweep = PETSC_FALSE; 10166dfbafefSToby Isaac fas->full_total = PETSC_FALSE; 101740244768SBarry Smith 101840244768SBarry Smith fas->eventsmoothsetup = 0; 101940244768SBarry Smith fas->eventsmoothsolve = 0; 102040244768SBarry Smith fas->eventresidual = 0; 102140244768SBarry Smith fas->eventinterprestrict = 0; 102240244768SBarry Smith PetscFunctionReturn(0); 102340244768SBarry Smith } 1024