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; 9f833ba53SLisandro Dalcin PetscErrorCode ierr; 10421d9b32SPeter Brune 11421d9b32SPeter Brune PetscFunctionBegin; 12ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr); 13ab8d36c9SPeter Brune ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr); 143dccd265SPeter Brune ierr = MatDestroy(&fas->inject);CHKERRQ(ierr); 15bccf9bb3SJed Brown ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr); 16bccf9bb3SJed Brown ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr); 17bccf9bb3SJed Brown ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr); 187cd33a7bSPeter Brune ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr); 197cd33a7bSPeter Brune ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr); 20d477d801SPeter Brune if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);} 21421d9b32SPeter Brune PetscFunctionReturn(0); 22421d9b32SPeter Brune } 23421d9b32SPeter Brune 2440244768SBarry Smith static PetscErrorCode SNESDestroy_FAS(SNES snes) 25421d9b32SPeter Brune { 26421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 27f833ba53SLisandro Dalcin PetscErrorCode ierr; 28421d9b32SPeter Brune 29421d9b32SPeter Brune PetscFunctionBegin; 30421d9b32SPeter Brune /* recursively resets and then destroys */ 31f833ba53SLisandro Dalcin ierr = SNESReset_FAS(snes);CHKERRQ(ierr); 321aa26658SKarl Rupp ierr = SNESDestroy(&fas->next);CHKERRQ(ierr); 33421d9b32SPeter Brune ierr = PetscFree(fas);CHKERRQ(ierr); 34421d9b32SPeter Brune PetscFunctionReturn(0); 35421d9b32SPeter Brune } 36421d9b32SPeter Brune 37f833ba53SLisandro Dalcin static PetscErrorCode SNESFASSetUpLineSearch_Private(SNES snes, SNES smooth) 38f833ba53SLisandro Dalcin { 39f833ba53SLisandro Dalcin SNESLineSearch linesearch; 40f833ba53SLisandro Dalcin SNESLineSearch slinesearch; 41f833ba53SLisandro Dalcin void *lsprectx,*lspostctx; 42f833ba53SLisandro Dalcin PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*); 43f833ba53SLisandro Dalcin PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*); 44f833ba53SLisandro Dalcin PetscErrorCode ierr; 45f833ba53SLisandro Dalcin 46f833ba53SLisandro Dalcin PetscFunctionBegin; 47f833ba53SLisandro Dalcin if (!snes->linesearch) PetscFunctionReturn(0); 48f833ba53SLisandro Dalcin ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr); 49f833ba53SLisandro Dalcin ierr = SNESGetLineSearch(smooth,&slinesearch);CHKERRQ(ierr); 50f833ba53SLisandro Dalcin ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr); 51f833ba53SLisandro Dalcin ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr); 52f833ba53SLisandro Dalcin ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr); 53f833ba53SLisandro Dalcin ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr); 54f833ba53SLisandro Dalcin ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr); 55f833ba53SLisandro Dalcin PetscFunctionReturn(0); 56f833ba53SLisandro Dalcin } 57f833ba53SLisandro Dalcin 58f833ba53SLisandro Dalcin static PetscErrorCode SNESFASCycleSetUpSmoother_Private(SNES snes, SNES smooth) 59f833ba53SLisandro Dalcin { 60f833ba53SLisandro Dalcin SNES_FAS *fas = (SNES_FAS*) snes->data; 61f833ba53SLisandro Dalcin PetscErrorCode ierr; 62f833ba53SLisandro Dalcin 63f833ba53SLisandro Dalcin PetscFunctionBegin; 64f833ba53SLisandro Dalcin ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)smooth);CHKERRQ(ierr); 65f833ba53SLisandro Dalcin ierr = SNESSetFromOptions(smooth);CHKERRQ(ierr); 66f833ba53SLisandro Dalcin ierr = SNESFASSetUpLineSearch_Private(snes, smooth);CHKERRQ(ierr); 67f833ba53SLisandro Dalcin 68f833ba53SLisandro Dalcin ierr = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr); 69f833ba53SLisandro Dalcin ierr = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr); 70f833ba53SLisandro Dalcin ierr = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr); 71f833ba53SLisandro Dalcin smooth->vec_sol = snes->vec_sol; 72f833ba53SLisandro Dalcin smooth->vec_sol_update = snes->vec_sol_update; 73f833ba53SLisandro Dalcin smooth->vec_func = snes->vec_func; 74f833ba53SLisandro Dalcin 75f833ba53SLisandro Dalcin if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,smooth,0,0,0);CHKERRQ(ierr);} 76f833ba53SLisandro Dalcin ierr = SNESSetUp(smooth);CHKERRQ(ierr); 77f833ba53SLisandro Dalcin if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,smooth,0,0,0);CHKERRQ(ierr);} 78f833ba53SLisandro Dalcin PetscFunctionReturn(0); 79f833ba53SLisandro Dalcin } 80f833ba53SLisandro Dalcin 8140244768SBarry Smith static PetscErrorCode SNESSetUp_FAS(SNES snes) 82421d9b32SPeter Brune { 8348bfdf8aSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 84421d9b32SPeter Brune PetscErrorCode ierr; 85d1adcc6fSPeter Brune PetscInt dm_levels; 86ab8d36c9SPeter Brune SNES next; 87f833ba53SLisandro Dalcin PetscBool isFine, hasCreateRestriction, hasCreateInjection; 88eff52c0eSPeter Brune 896b2b7091SBarry Smith PetscFunctionBegin; 90ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 91ab8d36c9SPeter Brune if (fas->usedmfornumberoflevels && isFine) { 92d1adcc6fSPeter Brune ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr); 93d1adcc6fSPeter Brune dm_levels++; 94cc05f883SPeter Brune if (dm_levels > fas->levels) { 953dccd265SPeter Brune /* reset the number of levels */ 960298fd71SBarry Smith ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr); 97cc05f883SPeter Brune ierr = SNESSetFromOptions(snes);CHKERRQ(ierr); 98d1adcc6fSPeter Brune } 99d1adcc6fSPeter Brune } 100ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 101ab8d36c9SPeter Brune if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */ 1023dccd265SPeter Brune 103fa0ddf94SBarry Smith ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */ 104cc05f883SPeter Brune 105ab8d36c9SPeter Brune /* set up the smoothers if they haven't already been set up */ 106ab8d36c9SPeter Brune if (!fas->smoothd) { 107ab8d36c9SPeter Brune ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr); 108ab8d36c9SPeter Brune } 109ab8d36c9SPeter Brune 11079d9a41aSPeter Brune if (snes->dm) { 111ab8d36c9SPeter Brune /* set the smoother DMs properly */ 112f833ba53SLisandro Dalcin if (fas->smoothu) {ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);} 113ab8d36c9SPeter Brune ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr); 11479d9a41aSPeter Brune /* construct EVERYTHING from the DM -- including the progressive set of smoothers */ 115ab8d36c9SPeter Brune if (next) { 11679d9a41aSPeter Brune /* for now -- assume the DM and the evaluation functions have been set externally */ 117ab8d36c9SPeter Brune if (!next->dm) { 118ce94432eSBarry Smith ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr); 119ab8d36c9SPeter Brune ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr); 12079d9a41aSPeter Brune } 12179d9a41aSPeter Brune /* set the interpolation and restriction from the DM */ 12279d9a41aSPeter Brune if (!fas->interpolate) { 123ab8d36c9SPeter Brune ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr); 124bccf9bb3SJed Brown if (!fas->restrct) { 1250a7266b2SPatrick Farrell ierr = DMHasCreateRestriction(next->dm, &hasCreateRestriction);CHKERRQ(ierr); 1260a7266b2SPatrick Farrell /* DM can create restrictions, use that */ 1270a7266b2SPatrick Farrell if (hasCreateRestriction) { 1280a7266b2SPatrick Farrell ierr = DMCreateRestriction(next->dm, snes->dm, &fas->restrct);CHKERRQ(ierr); 1290a7266b2SPatrick Farrell } else { 130bccf9bb3SJed Brown ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr); 13179d9a41aSPeter Brune fas->restrct = fas->interpolate; 13279d9a41aSPeter Brune } 133bccf9bb3SJed Brown } 1340a7266b2SPatrick Farrell } 13579d9a41aSPeter Brune /* set the injection from the DM */ 13679d9a41aSPeter Brune if (!fas->inject) { 137f833ba53SLisandro Dalcin ierr = DMHasCreateInjection(next->dm, &hasCreateInjection);CHKERRQ(ierr); 138f833ba53SLisandro Dalcin if (hasCreateInjection) { 1396dbf9973SLawrence Mitchell ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr); 14079d9a41aSPeter Brune } 14179d9a41aSPeter Brune } 14279d9a41aSPeter Brune } 14323e68893SLawrence Mitchell } 144f833ba53SLisandro Dalcin 14579d9a41aSPeter Brune /*pass the smoother, function, and jacobian up to the next level if it's not user set already */ 14679d9a41aSPeter Brune if (fas->galerkin) { 1471aa26658SKarl Rupp if (next) { 14825acbd8eSLisandro Dalcin ierr = SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next);CHKERRQ(ierr); 1491aa26658SKarl Rupp } 1501aa26658SKarl Rupp if (fas->smoothd && fas->level != fas->levels - 1) { 15125acbd8eSLisandro Dalcin ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes);CHKERRQ(ierr); 1521aa26658SKarl Rupp } 1531aa26658SKarl Rupp if (fas->smoothu && fas->level != fas->levels - 1) { 15425acbd8eSLisandro Dalcin ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes);CHKERRQ(ierr); 1551aa26658SKarl Rupp } 15679d9a41aSPeter Brune } 15779d9a41aSPeter Brune 158534ebe21SPeter Brune /* sets the down (pre) smoother's default norm and sets it from options */ 159534ebe21SPeter Brune if (fas->smoothd) { 160bc3f2f05SPeter Brune if (fas->level == 0 && fas->levels != 1) { 161365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr); 162534ebe21SPeter Brune } else { 163365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 164534ebe21SPeter Brune } 165f833ba53SLisandro Dalcin ierr = SNESFASCycleSetUpSmoother_Private(snes, fas->smoothd);CHKERRQ(ierr); 166534ebe21SPeter Brune } 167534ebe21SPeter Brune 168534ebe21SPeter Brune /* sets the up (post) smoother's default norm and sets it from options */ 169534ebe21SPeter Brune if (fas->smoothu) { 170534ebe21SPeter Brune if (fas->level != fas->levels - 1) { 171365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr); 172534ebe21SPeter Brune } else { 173365a6726SPeter Brune ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr); 174534ebe21SPeter Brune } 175f833ba53SLisandro Dalcin ierr = SNESFASCycleSetUpSmoother_Private(snes, fas->smoothu);CHKERRQ(ierr); 176534ebe21SPeter Brune } 177d06165b7SPeter Brune 178ab8d36c9SPeter Brune if (next) { 17979d9a41aSPeter Brune /* gotta set up the solution vector for this to work */ 180ab8d36c9SPeter Brune if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);} 181ab8d36c9SPeter Brune if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);} 1827fce8c19SPeter Brune ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr); 183f833ba53SLisandro Dalcin ierr = SNESFASSetUpLineSearch_Private(snes, next);CHKERRQ(ierr); 184ab8d36c9SPeter Brune ierr = SNESSetUp(next);CHKERRQ(ierr); 18579d9a41aSPeter Brune } 186f833ba53SLisandro Dalcin 1876273346dSPeter Brune /* setup FAS work vectors */ 1886273346dSPeter Brune if (fas->galerkin) { 1896273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr); 1906273346dSPeter Brune ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr); 1916273346dSPeter Brune } 192421d9b32SPeter Brune PetscFunctionReturn(0); 193421d9b32SPeter Brune } 194421d9b32SPeter Brune 19540244768SBarry Smith static PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes) 196421d9b32SPeter Brune { 197ee78dd50SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 198ee78dd50SPeter Brune PetscInt levels = 1; 19987f44e3fSPeter Brune PetscBool flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE; 200421d9b32SPeter Brune PetscErrorCode ierr; 20107144faaSPeter Brune SNESFASType fastype; 202fde0ff24SPeter Brune const char *optionsprefix; 203f1c6b773SPeter Brune SNESLineSearch linesearch; 20466585501SPeter Brune PetscInt m, n_up, n_down; 205ab8d36c9SPeter Brune SNES next; 206ab8d36c9SPeter Brune PetscBool isFine; 207421d9b32SPeter Brune 208421d9b32SPeter Brune PetscFunctionBegin; 209ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 210e55864a3SBarry Smith ierr = PetscOptionsHead(PetscOptionsObject,"SNESFAS Options-----------------------------------");CHKERRQ(ierr); 211ee78dd50SPeter Brune 212ab8d36c9SPeter Brune /* number of levels -- only process most options on the finest level */ 213ab8d36c9SPeter Brune if (isFine) { 214ee78dd50SPeter Brune ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr); 215c732cbdbSBarry Smith if (!flg && snes->dm) { 216c732cbdbSBarry Smith ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr); 217c732cbdbSBarry Smith levels++; 218d1adcc6fSPeter Brune fas->usedmfornumberoflevels = PETSC_TRUE; 219c732cbdbSBarry Smith } 2200298fd71SBarry Smith ierr = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr); 22107144faaSPeter Brune fastype = fas->fastype; 22207144faaSPeter Brune ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr); 22307144faaSPeter Brune if (flg) { 22407144faaSPeter Brune ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr); 22507144faaSPeter Brune } 226ee78dd50SPeter Brune 227fde0ff24SPeter Brune ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr); 228ab8d36c9SPeter Brune ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr); 229ab8d36c9SPeter Brune if (flg) { 230ab8d36c9SPeter Brune ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr); 231fde0ff24SPeter Brune } 23287f44e3fSPeter Brune ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr); 23387f44e3fSPeter Brune if (flg) { 23487f44e3fSPeter Brune ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr); 23587f44e3fSPeter Brune } 236fde0ff24SPeter Brune 237ab8d36c9SPeter Brune ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr); 238ab8d36c9SPeter Brune if (flg) { 239ab8d36c9SPeter Brune ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr); 240ab8d36c9SPeter Brune } 241ee78dd50SPeter Brune 242928e959bSPeter Brune if (fas->fastype == SNES_FAS_FULL) { 243ba672821SBarry Smith ierr = PetscOptionsBool("-snes_fas_full_downsweep","Smooth on the initial down sweep for full FAS cycles","SNESFASFullSetDownSweep",fas->full_downsweep,&fas->full_downsweep,&flg);CHKERRQ(ierr); 244ba672821SBarry Smith if (flg) {ierr = SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);} 245*6dfbafefSToby Isaac ierr = 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);CHKERRQ(ierr); 246*6dfbafefSToby Isaac if (flg) {ierr = SNESFASFullSetTotal(snes,fas->full_total);CHKERRQ(ierr);} 247928e959bSPeter Brune } 248928e959bSPeter Brune 24966585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr); 250162d76ddSPeter Brune 25166585501SPeter Brune ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr); 252162d76ddSPeter Brune 253d142ab34SLawrence Mitchell { 254d142ab34SLawrence Mitchell PetscViewer viewer; 255d142ab34SLawrence Mitchell PetscViewerFormat format; 25616413a6aSBarry Smith ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_fas_monitor",&viewer,&format,&monflg);CHKERRQ(ierr); 257d142ab34SLawrence Mitchell if (monflg) { 258d142ab34SLawrence Mitchell PetscViewerAndFormat *vf; 259d142ab34SLawrence Mitchell ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr); 260d142ab34SLawrence Mitchell ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr); 261d142ab34SLawrence Mitchell ierr = SNESFASSetMonitor(snes,vf,PETSC_TRUE);CHKERRQ(ierr); 262d142ab34SLawrence Mitchell } 263d142ab34SLawrence Mitchell } 2640dd27c6cSPeter Brune flg = PETSC_FALSE; 2650dd27c6cSPeter Brune monflg = PETSC_TRUE; 2660dd27c6cSPeter Brune ierr = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr); 2670dd27c6cSPeter Brune if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);} 268ab8d36c9SPeter Brune } 269ee78dd50SPeter Brune 270421d9b32SPeter Brune ierr = PetscOptionsTail();CHKERRQ(ierr); 271f833ba53SLisandro Dalcin 2728cc86e31SPeter Brune /* setup from the determined types if there is no pointwise procedure or smoother defined */ 273162d76ddSPeter Brune if (upflg) { 27466585501SPeter Brune ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr); 275162d76ddSPeter Brune } 276162d76ddSPeter Brune if (downflg) { 27766585501SPeter Brune ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr); 278162d76ddSPeter Brune } 279eff52c0eSPeter Brune 2809e764e56SPeter Brune /* set up the default line search for coarse grid corrections */ 2819e764e56SPeter Brune if (fas->fastype == SNES_FAS_ADDITIVE) { 2829e764e56SPeter Brune if (!snes->linesearch) { 2837601faf0SJed Brown ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr); 2841a4f838cSPeter Brune ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr); 2859e764e56SPeter Brune } 2869e764e56SPeter Brune } 2879e764e56SPeter Brune 288ee78dd50SPeter Brune /* recursive option setting for the smoothers */ 289f833ba53SLisandro Dalcin ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 290ab8d36c9SPeter Brune if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);} 291421d9b32SPeter Brune PetscFunctionReturn(0); 292421d9b32SPeter Brune } 293421d9b32SPeter Brune 2949804daf3SBarry Smith #include <petscdraw.h> 29540244768SBarry Smith static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer) 296421d9b32SPeter Brune { 297421d9b32SPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 298656ede7eSPeter Brune PetscBool isFine,iascii,isdraw; 299ab8d36c9SPeter Brune PetscInt i; 300421d9b32SPeter Brune PetscErrorCode ierr; 301ab8d36c9SPeter Brune SNES smoothu, smoothd, levelsnes; 302421d9b32SPeter Brune 303421d9b32SPeter Brune PetscFunctionBegin; 304ab8d36c9SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 305ab8d36c9SPeter Brune if (isFine) { 306251f4c67SDmitry Karpeev ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 307656ede7eSPeter Brune ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr); 308421d9b32SPeter Brune if (iascii) { 309efd4aadfSBarry Smith ierr = PetscViewerASCIIPrintf(viewer, " type is %s, levels=%D, cycles=%D\n", SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr); 310ab8d36c9SPeter Brune if (fas->galerkin) { 311ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 312421d9b32SPeter Brune } else { 313ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr); 314421d9b32SPeter Brune } 315ab8d36c9SPeter Brune for (i=0; i<fas->levels; i++) { 316ab8d36c9SPeter Brune ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr); 317ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr); 318ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr); 319ab8d36c9SPeter Brune if (!i) { 320ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr); 321421d9b32SPeter Brune } else { 322ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 323421d9b32SPeter Brune } 324ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 325166b3ea4SJed Brown if (smoothd) { 326ab8d36c9SPeter Brune ierr = SNESView(smoothd,viewer);CHKERRQ(ierr); 327166b3ea4SJed Brown } else { 328166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 329166b3ea4SJed Brown } 330ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 331ab8d36c9SPeter Brune if (i && (smoothd == smoothu)) { 332ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr); 333ab8d36c9SPeter Brune } else if (i) { 334ab8d36c9SPeter Brune ierr = PetscViewerASCIIPrintf(viewer," Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr); 335ab8d36c9SPeter Brune ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr); 336166b3ea4SJed Brown if (smoothu) { 337ab8d36c9SPeter Brune ierr = SNESView(smoothu,viewer);CHKERRQ(ierr); 338166b3ea4SJed Brown } else { 339166b3ea4SJed Brown ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr); 340166b3ea4SJed Brown } 341ab8d36c9SPeter Brune ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr); 342ab8d36c9SPeter Brune } 343ab8d36c9SPeter Brune } 344656ede7eSPeter Brune } else if (isdraw) { 345656ede7eSPeter Brune PetscDraw draw; 346b4375e8dSPeter Brune PetscReal x,w,y,bottom,th,wth; 347656ede7eSPeter Brune SNES_FAS *curfas = fas; 348656ede7eSPeter Brune ierr = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr); 349656ede7eSPeter Brune ierr = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr); 350656ede7eSPeter Brune ierr = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr); 351656ede7eSPeter Brune bottom = y - th; 352656ede7eSPeter Brune while (curfas) { 353b4375e8dSPeter Brune if (!curfas->smoothu) { 354656ede7eSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr); 355656ede7eSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 356656ede7eSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 357b4375e8dSPeter Brune } else { 358b4375e8dSPeter Brune w = 0.5*PetscMin(1.0-x,x); 359b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr); 360b4375e8dSPeter Brune if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr); 361b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 362b4375e8dSPeter Brune ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr); 363b4375e8dSPeter Brune if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr); 364b4375e8dSPeter Brune ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr); 365b4375e8dSPeter Brune } 366656ede7eSPeter Brune /* this is totally bogus but we have no way of knowing how low the previous one was draw to */ 367656ede7eSPeter Brune bottom -= 5*th; 3681aa26658SKarl Rupp if (curfas->next) curfas = (SNES_FAS*)curfas->next->data; 3690298fd71SBarry Smith else curfas = NULL; 370656ede7eSPeter Brune } 371421d9b32SPeter Brune } 372ab8d36c9SPeter Brune } 373421d9b32SPeter Brune PetscFunctionReturn(0); 374421d9b32SPeter Brune } 375421d9b32SPeter Brune 37639bd7f45SPeter Brune /* 37739bd7f45SPeter Brune Defines the action of the downsmoother 37839bd7f45SPeter Brune */ 37940244768SBarry Smith static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 380b9c2fdf1SPeter Brune { 381f833ba53SLisandro Dalcin PetscErrorCode ierr; 382742fe5e2SPeter Brune SNESConvergedReason reason; 383ab8d36c9SPeter Brune Vec FPC; 384ab8d36c9SPeter Brune SNES smoothd; 3856cbb2f26SLawrence Mitchell PetscBool flg; 3860dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 3876e111a19SKarl Rupp 388421d9b32SPeter Brune PetscFunctionBegin; 389ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr); 390e4ed7901SPeter Brune ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr); 391217044c2SLisandro Dalcin if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,smoothd,B,X,0);CHKERRQ(ierr);} 392ab8d36c9SPeter Brune ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr); 393217044c2SLisandro Dalcin if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,smoothd,B,X,0);CHKERRQ(ierr);} 394742fe5e2SPeter Brune /* check convergence reason for the smoother */ 395ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr); 3961ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 397742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 398742fe5e2SPeter Brune PetscFunctionReturn(0); 399742fe5e2SPeter Brune } 4006cbb2f26SLawrence Mitchell 4010298fd71SBarry Smith ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr); 4026cbb2f26SLawrence Mitchell ierr = SNESGetAlwaysComputesFinalResidual(smoothd, &flg);CHKERRQ(ierr); 4036cbb2f26SLawrence Mitchell if (!flg) { 4046cbb2f26SLawrence Mitchell ierr = SNESComputeFunction(smoothd, X, FPC);CHKERRQ(ierr); 4056cbb2f26SLawrence Mitchell } 4064b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 40771dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 40839bd7f45SPeter Brune PetscFunctionReturn(0); 40939bd7f45SPeter Brune } 41039bd7f45SPeter Brune 41139bd7f45SPeter Brune 41239bd7f45SPeter Brune /* 41307144faaSPeter Brune Defines the action of the upsmoother 41439bd7f45SPeter Brune */ 41540244768SBarry Smith static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) 4160adebc6cSBarry Smith { 417f833ba53SLisandro Dalcin PetscErrorCode ierr; 41839bd7f45SPeter Brune SNESConvergedReason reason; 419ab8d36c9SPeter Brune Vec FPC; 420ab8d36c9SPeter Brune SNES smoothu; 4216cbb2f26SLawrence Mitchell PetscBool flg; 4220dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*) snes->data; 423ab8d36c9SPeter Brune 4246e111a19SKarl Rupp PetscFunctionBegin; 425ab8d36c9SPeter Brune ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr); 426217044c2SLisandro Dalcin if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,smoothu,0,0,0);CHKERRQ(ierr);} 427ab8d36c9SPeter Brune ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr); 428217044c2SLisandro Dalcin if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,smoothu,0,0,0);CHKERRQ(ierr);} 42939bd7f45SPeter Brune /* check convergence reason for the smoother */ 430ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr); 4311ff805c3SPeter Brune if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) { 43239bd7f45SPeter Brune snes->reason = SNES_DIVERGED_INNER; 43339bd7f45SPeter Brune PetscFunctionReturn(0); 43439bd7f45SPeter Brune } 4350298fd71SBarry Smith ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr); 4366cbb2f26SLawrence Mitchell ierr = SNESGetAlwaysComputesFinalResidual(smoothu, &flg);CHKERRQ(ierr); 4376cbb2f26SLawrence Mitchell if (!flg) { 4386cbb2f26SLawrence Mitchell ierr = SNESComputeFunction(smoothu, X, FPC);CHKERRQ(ierr); 4396cbb2f26SLawrence Mitchell } 4404b32a720SPeter Brune ierr = VecCopy(FPC, F);CHKERRQ(ierr); 44171dbe336SPeter Brune if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);} 44239bd7f45SPeter Brune PetscFunctionReturn(0); 44339bd7f45SPeter Brune } 44439bd7f45SPeter Brune 445938e4a01SJed Brown /*@ 446938e4a01SJed Brown SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level 447938e4a01SJed Brown 448938e4a01SJed Brown Collective 449938e4a01SJed Brown 450938e4a01SJed Brown Input Arguments: 451938e4a01SJed Brown . snes - SNESFAS 452938e4a01SJed Brown 453938e4a01SJed Brown Output Arguments: 454938e4a01SJed Brown . Xcoarse - vector on level one coarser than snes 455938e4a01SJed Brown 456938e4a01SJed Brown Level: developer 457938e4a01SJed Brown 458938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict() 459938e4a01SJed Brown @*/ 460938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse) 461938e4a01SJed Brown { 462938e4a01SJed Brown PetscErrorCode ierr; 463f833ba53SLisandro Dalcin SNES_FAS *fas; 464938e4a01SJed Brown 465938e4a01SJed Brown PetscFunctionBegin; 466f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS); 467f833ba53SLisandro Dalcin PetscValidPointer(Xcoarse,2); 468f833ba53SLisandro Dalcin fas = (SNES_FAS*)snes->data; 4691aa26658SKarl Rupp if (fas->rscale) { 4701aa26658SKarl Rupp ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr); 471f5af7f23SKarl Rupp } else if (fas->inject) { 4722a7a6963SBarry Smith ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr); 47313903a91SSatish Balay } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection"); 474938e4a01SJed Brown PetscFunctionReturn(0); 475938e4a01SJed Brown } 476938e4a01SJed Brown 477e9923e8dSJed Brown /*@ 478e9923e8dSJed Brown SNESFASRestrict - restrict a Vec to the next coarser level 479e9923e8dSJed Brown 480e9923e8dSJed Brown Collective 481e9923e8dSJed Brown 482e9923e8dSJed Brown Input Arguments: 483e9923e8dSJed Brown + fine - SNES from which to restrict 484e9923e8dSJed Brown - Xfine - vector to restrict 485e9923e8dSJed Brown 486e9923e8dSJed Brown Output Arguments: 487e9923e8dSJed Brown . Xcoarse - result of restriction 488e9923e8dSJed Brown 489e9923e8dSJed Brown Level: developer 490e9923e8dSJed Brown 491e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection() 492e9923e8dSJed Brown @*/ 493e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse) 494e9923e8dSJed Brown { 495e9923e8dSJed Brown PetscErrorCode ierr; 496f833ba53SLisandro Dalcin SNES_FAS *fas; 497e9923e8dSJed Brown 498e9923e8dSJed Brown PetscFunctionBegin; 499f833ba53SLisandro Dalcin PetscValidHeaderSpecificType(fine,SNES_CLASSID,1,SNESFAS); 500e9923e8dSJed Brown PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2); 501e9923e8dSJed Brown PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3); 502f833ba53SLisandro Dalcin fas = (SNES_FAS*)fine->data; 503e9923e8dSJed Brown if (fas->inject) { 504e9923e8dSJed Brown ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr); 505e9923e8dSJed Brown } else { 506e9923e8dSJed Brown ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr); 507e9923e8dSJed Brown ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr); 508e9923e8dSJed Brown } 509e9923e8dSJed Brown PetscFunctionReturn(0); 510e9923e8dSJed Brown } 511e9923e8dSJed Brown 51239bd7f45SPeter Brune /* 51339bd7f45SPeter Brune 514*6dfbafefSToby Isaac Performs a variant of FAS using the interpolated total coarse solution 515*6dfbafefSToby Isaac 516*6dfbafefSToby Isaac fine problem: F(x) = b 517*6dfbafefSToby Isaac coarse problem: F^c(x^c) = Rb, Initial guess Rx 518*6dfbafefSToby Isaac interpolated solution: x^f = I x^c (total solution interpolation 519*6dfbafefSToby Isaac 520*6dfbafefSToby Isaac */ 521*6dfbafefSToby Isaac static PetscErrorCode SNESFASInterpolatedCoarseSolution(SNES snes, Vec X, Vec X_new) 522*6dfbafefSToby Isaac { 523*6dfbafefSToby Isaac PetscErrorCode ierr; 524*6dfbafefSToby Isaac Vec X_c, B_c; 525*6dfbafefSToby Isaac SNESConvergedReason reason; 526*6dfbafefSToby Isaac SNES next; 527*6dfbafefSToby Isaac Mat restrct, interpolate; 528*6dfbafefSToby Isaac SNES_FAS *fasc; 529*6dfbafefSToby Isaac 530*6dfbafefSToby Isaac PetscFunctionBegin; 531*6dfbafefSToby Isaac ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 532*6dfbafefSToby Isaac if (next) { 533*6dfbafefSToby Isaac fasc = (SNES_FAS*)next->data; 534*6dfbafefSToby Isaac 535*6dfbafefSToby Isaac ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 536*6dfbafefSToby Isaac ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 537*6dfbafefSToby Isaac 538*6dfbafefSToby Isaac X_c = next->vec_sol; 539*6dfbafefSToby Isaac 540*6dfbafefSToby Isaac if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);} 541*6dfbafefSToby Isaac /* restrict the total solution: Rb */ 542*6dfbafefSToby Isaac ierr = SNESFASRestrict(snes, X, X_c);CHKERRQ(ierr); 543*6dfbafefSToby Isaac B_c = next->vec_rhs; 544*6dfbafefSToby Isaac if (snes->vec_rhs) { 545*6dfbafefSToby Isaac /* restrict the total rhs defect: Rb */ 546*6dfbafefSToby Isaac ierr = MatRestrict(restrct, snes->vec_rhs, B_c);CHKERRQ(ierr); 547*6dfbafefSToby Isaac } else { 548*6dfbafefSToby Isaac ierr = VecSet(B_c, 0.);CHKERRQ(ierr); 549*6dfbafefSToby Isaac } 550*6dfbafefSToby Isaac if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);} 551*6dfbafefSToby Isaac 552*6dfbafefSToby Isaac ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 553*6dfbafefSToby Isaac ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 554*6dfbafefSToby Isaac if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 555*6dfbafefSToby Isaac snes->reason = SNES_DIVERGED_INNER; 556*6dfbafefSToby Isaac PetscFunctionReturn(0); 557*6dfbafefSToby Isaac } 558*6dfbafefSToby Isaac /* x^f <- Ix^c*/ 559*6dfbafefSToby Isaac DM dmc,dmf; 560*6dfbafefSToby Isaac 561*6dfbafefSToby Isaac ierr = SNESGetDM(next, &dmc);CHKERRQ(ierr); 562*6dfbafefSToby Isaac ierr = SNESGetDM(snes, &dmf);CHKERRQ(ierr); 563*6dfbafefSToby Isaac if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);} 564*6dfbafefSToby Isaac ierr = DMInterpolateSolution(dmc, dmf, interpolate, X_c, X_new);CHKERRQ(ierr); 565*6dfbafefSToby Isaac if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);} 566*6dfbafefSToby Isaac ierr = PetscObjectSetName((PetscObject) X_c, "Coarse solution");CHKERRQ(ierr); 567*6dfbafefSToby Isaac ierr = VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view");CHKERRQ(ierr); 568*6dfbafefSToby Isaac ierr = PetscObjectSetName((PetscObject) X_new, "Updated Fine solution");CHKERRQ(ierr); 569*6dfbafefSToby Isaac ierr = VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view");CHKERRQ(ierr); 570*6dfbafefSToby Isaac } 571*6dfbafefSToby Isaac PetscFunctionReturn(0); 572*6dfbafefSToby Isaac } 573*6dfbafefSToby Isaac 574*6dfbafefSToby Isaac /* 575*6dfbafefSToby Isaac 57639bd7f45SPeter Brune Performs the FAS coarse correction as: 57739bd7f45SPeter Brune 578b5527d98SMatthew G. Knepley fine problem: F(x) = b 579b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c 58039bd7f45SPeter Brune 581b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 58239bd7f45SPeter Brune 58339bd7f45SPeter Brune */ 5840adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) 5850adebc6cSBarry Smith { 58639bd7f45SPeter Brune PetscErrorCode ierr; 58739bd7f45SPeter Brune Vec X_c, Xo_c, F_c, B_c; 58839bd7f45SPeter Brune SNESConvergedReason reason; 589ab8d36c9SPeter Brune SNES next; 590ab8d36c9SPeter Brune Mat restrct, interpolate; 5910dd27c6cSPeter Brune SNES_FAS *fasc; 5925fd66863SKarl Rupp 59339bd7f45SPeter Brune PetscFunctionBegin; 594ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 595ab8d36c9SPeter Brune if (next) { 5960dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 5970dd27c6cSPeter Brune 598ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 599ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 600ab8d36c9SPeter Brune 601ab8d36c9SPeter Brune X_c = next->vec_sol; 602ab8d36c9SPeter Brune Xo_c = next->work[0]; 603ab8d36c9SPeter Brune F_c = next->vec_func; 604ab8d36c9SPeter Brune B_c = next->vec_rhs; 605efe1f98aSPeter Brune 606217044c2SLisandro Dalcin if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);} 607938e4a01SJed Brown ierr = SNESFASRestrict(snes, X, Xo_c);CHKERRQ(ierr); 6085609cbf2SMatthew G. Knepley /* restrict the defect: R(F(x) - b) */ 609ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 610217044c2SLisandro Dalcin if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);} 6110dd27c6cSPeter Brune 612217044c2SLisandro Dalcin if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);} 6135609cbf2SMatthew G. Knepley /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */ 614ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 615217044c2SLisandro Dalcin if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);} 6160dd27c6cSPeter Brune 6170dd27c6cSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */ 618e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 619b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 620e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 621ee78dd50SPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 622ee78dd50SPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 623c90fad12SPeter Brune 624c90fad12SPeter Brune /* recurse to the next level */ 625e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 626ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 627ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 628742fe5e2SPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 629742fe5e2SPeter Brune snes->reason = SNES_DIVERGED_INNER; 630742fe5e2SPeter Brune PetscFunctionReturn(0); 631742fe5e2SPeter Brune } 632fa9694d7SPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 633fa9694d7SPeter Brune ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr); 6340dd27c6cSPeter Brune 635217044c2SLisandro Dalcin if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);} 636ab8d36c9SPeter Brune ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr); 637*6dfbafefSToby Isaac if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);} 638fb8dc43dSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) X_c, "Coarse correction");CHKERRQ(ierr); 639fb8dc43dSMatthew G. Knepley ierr = VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view");CHKERRQ(ierr); 640fb8dc43dSMatthew G. Knepley ierr = PetscObjectSetName((PetscObject) X_new, "Updated Fine solution");CHKERRQ(ierr); 641fb8dc43dSMatthew G. Knepley ierr = VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view");CHKERRQ(ierr); 642293a7e31SPeter Brune } 64339bd7f45SPeter Brune PetscFunctionReturn(0); 64439bd7f45SPeter Brune } 64539bd7f45SPeter Brune 64639bd7f45SPeter Brune /* 64739bd7f45SPeter Brune 64839bd7f45SPeter Brune The additive cycle looks like: 64939bd7f45SPeter Brune 65007144faaSPeter Brune xhat = x 65107144faaSPeter Brune xhat = dS(x, b) 65207144faaSPeter Brune x = coarsecorrection(xhat, b_d) 65307144faaSPeter Brune x = x + nu*(xhat - x); 65439bd7f45SPeter Brune (optional) x = uS(x, b) 65539bd7f45SPeter Brune 65639bd7f45SPeter Brune With the coarse RHS (defect correction) as below. 65739bd7f45SPeter Brune 65839bd7f45SPeter Brune */ 65940244768SBarry Smith static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) 6600adebc6cSBarry Smith { 66107144faaSPeter Brune Vec F, B, Xhat; 66222c1e704SPeter Brune Vec X_c, Xo_c, F_c, B_c; 66339bd7f45SPeter Brune PetscErrorCode ierr; 66407144faaSPeter Brune SNESConvergedReason reason; 66522c1e704SPeter Brune PetscReal xnorm, fnorm, ynorm; 666422a814eSBarry Smith SNESLineSearchReason lsresult; 667ab8d36c9SPeter Brune SNES next; 668ab8d36c9SPeter Brune Mat restrct, interpolate; 6690dd27c6cSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data,*fasc; 6700dd27c6cSPeter Brune 67139bd7f45SPeter Brune PetscFunctionBegin; 672ab8d36c9SPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 67339bd7f45SPeter Brune F = snes->vec_func; 67439bd7f45SPeter Brune B = snes->vec_rhs; 675e7f468e7SPeter Brune Xhat = snes->work[1]; 67607144faaSPeter Brune ierr = VecCopy(X, Xhat);CHKERRQ(ierr); 67707144faaSPeter Brune /* recurse first */ 678ab8d36c9SPeter Brune if (next) { 6790dd27c6cSPeter Brune fasc = (SNES_FAS*)next->data; 680ab8d36c9SPeter Brune ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr); 681ab8d36c9SPeter Brune ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr); 682217044c2SLisandro Dalcin if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);} 68307144faaSPeter Brune ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr); 684217044c2SLisandro Dalcin if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);} 685c2a02606SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); 686ab8d36c9SPeter Brune X_c = next->vec_sol; 687ab8d36c9SPeter Brune Xo_c = next->work[0]; 688ab8d36c9SPeter Brune F_c = next->vec_func; 689ab8d36c9SPeter Brune B_c = next->vec_rhs; 69039bd7f45SPeter Brune 691938e4a01SJed Brown ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr); 69207144faaSPeter Brune /* restrict the defect */ 693ab8d36c9SPeter Brune ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr); 69407144faaSPeter Brune 69507144faaSPeter Brune /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */ 696217044c2SLisandro Dalcin if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);} 697ab8d36c9SPeter Brune ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr); 698217044c2SLisandro Dalcin if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);} 699e4ed7901SPeter Brune ierr = VecCopy(B_c, X_c);CHKERRQ(ierr); 700b9c2fdf1SPeter Brune ierr = VecCopy(F_c, B_c);CHKERRQ(ierr); 701e4ed7901SPeter Brune ierr = VecCopy(X_c, F_c);CHKERRQ(ierr); 70207144faaSPeter Brune /* set initial guess of the coarse problem to the projected fine solution */ 70307144faaSPeter Brune ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr); 70407144faaSPeter Brune 70507144faaSPeter Brune /* recurse */ 706e4ed7901SPeter Brune ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr); 707ab8d36c9SPeter Brune ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr); 70807144faaSPeter Brune 70907144faaSPeter Brune /* smooth on this level */ 71091f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr); 71107144faaSPeter Brune 712ab8d36c9SPeter Brune ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr); 71307144faaSPeter Brune if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) { 71407144faaSPeter Brune snes->reason = SNES_DIVERGED_INNER; 71507144faaSPeter Brune PetscFunctionReturn(0); 71607144faaSPeter Brune } 71707144faaSPeter Brune 71807144faaSPeter Brune /* correct as x <- x + I(x^c - Rx)*/ 719c68acad4SPeter Brune ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr); 720ab8d36c9SPeter Brune ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr); 72107144faaSPeter Brune 722ddebd997SPeter Brune /* additive correction of the coarse direction*/ 723f1c6b773SPeter Brune ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr); 724422a814eSBarry Smith ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr); 725422a814eSBarry Smith ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr); 726422a814eSBarry Smith if (lsresult) { 7279e764e56SPeter Brune if (++snes->numFailures >= snes->maxFailures) { 7289e764e56SPeter Brune snes->reason = SNES_DIVERGED_LINE_SEARCH; 7299e764e56SPeter Brune PetscFunctionReturn(0); 7309e764e56SPeter Brune } 7319e764e56SPeter Brune } 73207144faaSPeter Brune } else { 73391f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 73407144faaSPeter Brune } 73539bd7f45SPeter Brune PetscFunctionReturn(0); 73639bd7f45SPeter Brune } 73739bd7f45SPeter Brune 73839bd7f45SPeter Brune /* 73939bd7f45SPeter Brune 74039bd7f45SPeter Brune Defines the FAS cycle as: 74139bd7f45SPeter Brune 742b5527d98SMatthew G. Knepley fine problem: F(x) = b 74339bd7f45SPeter Brune coarse problem: F^c(x) = b^c 74439bd7f45SPeter Brune 745b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b) 74639bd7f45SPeter Brune 74739bd7f45SPeter Brune correction: 74839bd7f45SPeter Brune 74939bd7f45SPeter Brune x = x + I(x^c - Rx) 75039bd7f45SPeter Brune 75139bd7f45SPeter Brune */ 75240244768SBarry Smith static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) 7530adebc6cSBarry Smith { 75439bd7f45SPeter Brune 75539bd7f45SPeter Brune PetscErrorCode ierr; 75639bd7f45SPeter Brune Vec F,B; 75734d65b3cSPeter Brune SNES next; 75839bd7f45SPeter Brune 75939bd7f45SPeter Brune PetscFunctionBegin; 76039bd7f45SPeter Brune F = snes->vec_func; 76139bd7f45SPeter Brune B = snes->vec_rhs; 76239bd7f45SPeter Brune /* pre-smooth -- just update using the pre-smoother */ 76334d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr); 76491f99d7cSPeter Brune ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 76534d65b3cSPeter Brune if (next) { 7668c40d5fbSBarry Smith ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr); 76791f99d7cSPeter Brune ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr); 768fe6f9142SPeter Brune } 769fa9694d7SPeter Brune PetscFunctionReturn(0); 770421d9b32SPeter Brune } 771421d9b32SPeter Brune 77240244768SBarry Smith static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes) 7738c94862eSPeter Brune { 7748c94862eSPeter Brune SNES next; 7758c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 7768c94862eSPeter Brune PetscBool isFine; 7778c94862eSPeter Brune PetscErrorCode ierr; 7788c94862eSPeter Brune 7798c94862eSPeter Brune PetscFunctionBegin; 7808c94862eSPeter Brune /* pre-smooth -- just update using the pre-smoother */ 7818c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 7828c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 7838c94862eSPeter Brune fas->full_stage = 0; 7848c94862eSPeter Brune if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);} 7858c94862eSPeter Brune PetscFunctionReturn(0); 7868c94862eSPeter Brune } 7878c94862eSPeter Brune 78840244768SBarry Smith static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X) 7898c94862eSPeter Brune { 7908c94862eSPeter Brune PetscErrorCode ierr; 7918c94862eSPeter Brune Vec F,B; 7928c94862eSPeter Brune SNES_FAS *fas = (SNES_FAS*)snes->data; 7938c94862eSPeter Brune PetscBool isFine; 7948c94862eSPeter Brune SNES next; 7958c94862eSPeter Brune 7968c94862eSPeter Brune PetscFunctionBegin; 7978c94862eSPeter Brune F = snes->vec_func; 7988c94862eSPeter Brune B = snes->vec_rhs; 7998c94862eSPeter Brune ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr); 8008c94862eSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 8018c94862eSPeter Brune 8028c94862eSPeter Brune if (isFine) { 8038c94862eSPeter Brune ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr); 8048c94862eSPeter Brune } 8058c94862eSPeter Brune 8068c94862eSPeter Brune if (fas->full_stage == 0) { 807928e959bSPeter Brune /* downsweep */ 8088c94862eSPeter Brune if (next) { 8098c94862eSPeter Brune if (fas->level != 1) next->max_its += 1; 810e140b65aSPatrick Farrell if (fas->full_downsweep) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 811e140b65aSPatrick Farrell fas->full_downsweep = PETSC_TRUE; 812*6dfbafefSToby Isaac if (fas->full_total) {ierr = SNESFASInterpolatedCoarseSolution(snes,X,X);CHKERRQ(ierr);} 813*6dfbafefSToby Isaac else {ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);} 814*6dfbafefSToby Isaac fas->full_total = PETSC_FALSE; 8158c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8168c94862eSPeter Brune if (fas->level != 1) next->max_its -= 1; 8178c94862eSPeter Brune } else { 818a3a80b83SMatthew G. Knepley /* The smoother on the coarse level is the coarse solver */ 8198c94862eSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8208c94862eSPeter Brune } 8218c94862eSPeter Brune fas->full_stage = 1; 8228c94862eSPeter Brune } else if (fas->full_stage == 1) { 8238c94862eSPeter Brune if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);} 8248c94862eSPeter Brune if (next) { 8258c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8268c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8278c94862eSPeter Brune } 8288c94862eSPeter Brune } 8298c94862eSPeter Brune /* final v-cycle */ 8308c94862eSPeter Brune if (isFine) { 8318c94862eSPeter Brune if (next) { 8328c94862eSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 8338c94862eSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 8348c94862eSPeter Brune } 8358c94862eSPeter Brune } 8368c94862eSPeter Brune PetscFunctionReturn(0); 8378c94862eSPeter Brune } 8388c94862eSPeter Brune 83940244768SBarry Smith static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X) 84034d65b3cSPeter Brune { 84134d65b3cSPeter Brune PetscErrorCode ierr; 84234d65b3cSPeter Brune Vec F,B; 84334d65b3cSPeter Brune SNES next; 84434d65b3cSPeter Brune 84534d65b3cSPeter Brune PetscFunctionBegin; 84634d65b3cSPeter Brune F = snes->vec_func; 84734d65b3cSPeter Brune B = snes->vec_rhs; 84834d65b3cSPeter Brune ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr); 84934d65b3cSPeter Brune if (next) { 85034d65b3cSPeter Brune ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr); 85134d65b3cSPeter Brune ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 85234d65b3cSPeter Brune } else { 85334d65b3cSPeter Brune ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr); 85434d65b3cSPeter Brune } 85534d65b3cSPeter Brune PetscFunctionReturn(0); 85634d65b3cSPeter Brune } 85734d65b3cSPeter Brune 858fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE; 859fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n" 860fffbeea8SBarry Smith " title = {Composing Scalable Nonlinear Algebraic Solvers},\n" 861fffbeea8SBarry Smith " author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n" 862fffbeea8SBarry Smith " year = 2013,\n" 863fffbeea8SBarry Smith " type = Preprint,\n" 864fffbeea8SBarry Smith " number = {ANL/MCS-P2010-0112},\n" 865fffbeea8SBarry Smith " institution = {Argonne National Laboratory}\n}\n"; 866fffbeea8SBarry Smith 86740244768SBarry Smith static PetscErrorCode SNESSolve_FAS(SNES snes) 868421d9b32SPeter Brune { 869fa9694d7SPeter Brune PetscErrorCode ierr; 870f833ba53SLisandro Dalcin PetscInt i; 871ddb5aff1SPeter Brune Vec X, F; 872fe6f9142SPeter Brune PetscReal fnorm; 873b17ce1afSJed Brown SNES_FAS *fas = (SNES_FAS*)snes->data,*ffas; 874b17ce1afSJed Brown DM dm; 875e70c42e5SPeter Brune PetscBool isFine; 876b17ce1afSJed Brown 877421d9b32SPeter Brune PetscFunctionBegin; 8786c4ed002SBarry Smith if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name); 879c579b300SPatrick Farrell 880fffbeea8SBarry Smith ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr); 881fe6f9142SPeter Brune snes->reason = SNES_CONVERGED_ITERATING; 882fa9694d7SPeter Brune X = snes->vec_sol; 883f5a6d4f9SBarry Smith F = snes->vec_func; 884293a7e31SPeter Brune 88518a66777SPeter Brune ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr); 886293a7e31SPeter Brune /* norm setup */ 887e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 888fe6f9142SPeter Brune snes->iter = 0; 889f833ba53SLisandro Dalcin snes->norm = 0; 890e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 891e4ed7901SPeter Brune if (!snes->vec_func_init_set) { 892217044c2SLisandro Dalcin if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);} 893fe6f9142SPeter Brune ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr); 894217044c2SLisandro Dalcin if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);} 8951aa26658SKarl Rupp } else snes->vec_func_init_set = PETSC_FALSE; 896e4ed7901SPeter Brune 897fe6f9142SPeter Brune ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F|| */ 898422a814eSBarry Smith SNESCheckFunctionNorm(snes,fnorm); 899e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 900fe6f9142SPeter Brune snes->norm = fnorm; 901e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 902a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr); 903f833ba53SLisandro Dalcin ierr = SNESMonitor(snes,snes->iter,fnorm);CHKERRQ(ierr); 904fe6f9142SPeter Brune 905fe6f9142SPeter Brune /* test convergence */ 906fe6f9142SPeter Brune ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 907fe6f9142SPeter Brune if (snes->reason) PetscFunctionReturn(0); 908e4ed7901SPeter Brune 909b17ce1afSJed Brown 910b9c2fdf1SPeter Brune if (isFine) { 911b9c2fdf1SPeter Brune /* propagate scale-dependent data up the hierarchy */ 912b17ce1afSJed Brown ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr); 913b17ce1afSJed Brown for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) { 914b17ce1afSJed Brown DM dmcoarse; 915b17ce1afSJed Brown ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr); 916b17ce1afSJed Brown ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr); 917b17ce1afSJed Brown dm = dmcoarse; 918b17ce1afSJed Brown } 919b9c2fdf1SPeter Brune } 920b17ce1afSJed Brown 921f833ba53SLisandro Dalcin for (i = 0; i < snes->max_its; i++) { 922fe6f9142SPeter Brune /* Call general purpose update function */ 923fe6f9142SPeter Brune if (snes->ops->update) { 924fe6f9142SPeter Brune ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr); 925fe6f9142SPeter Brune } 926f833ba53SLisandro Dalcin 92707144faaSPeter Brune if (fas->fastype == SNES_FAS_MULTIPLICATIVE) { 92891f99d7cSPeter Brune ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr); 9298c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_ADDITIVE) { 93091f99d7cSPeter Brune ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr); 9318c94862eSPeter Brune } else if (fas->fastype == SNES_FAS_FULL) { 9328c94862eSPeter Brune ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr); 93334d65b3cSPeter Brune } else if (fas->fastype == SNES_FAS_KASKADE) { 93434d65b3cSPeter Brune ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr); 9356c4ed002SBarry Smith } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type"); 936742fe5e2SPeter Brune 937742fe5e2SPeter Brune /* check for FAS cycle divergence */ 9381aa26658SKarl Rupp if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0); 939b9c2fdf1SPeter Brune 940c90fad12SPeter Brune /* Monitor convergence */ 941e04113cfSBarry Smith ierr = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr); 942c90fad12SPeter Brune snes->iter = i+1; 943e04113cfSBarry Smith ierr = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr); 944a71f0d7dSBarry Smith ierr = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr); 945c90fad12SPeter Brune ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr); 946c90fad12SPeter Brune /* Test for convergence */ 94766585501SPeter Brune if (isFine) { 948b9c2fdf1SPeter Brune ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr); 949c90fad12SPeter Brune if (snes->reason) break; 950fe6f9142SPeter Brune } 95166585501SPeter Brune } 952f833ba53SLisandro Dalcin if (i == snes->max_its) { 953f833ba53SLisandro Dalcin ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", i);CHKERRQ(ierr); 954fe6f9142SPeter Brune if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT; 955fe6f9142SPeter Brune } 956421d9b32SPeter Brune PetscFunctionReturn(0); 957421d9b32SPeter Brune } 95840244768SBarry Smith 95940244768SBarry Smith /*MC 96040244768SBarry Smith 96140244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver. 96240244768SBarry Smith 96340244768SBarry Smith The nonlinear problem is solved by correction using coarse versions 96440244768SBarry Smith of the nonlinear problem. This problem is perturbed so that a projected 96540244768SBarry Smith solution of the fine problem elicits no correction from the coarse problem. 96640244768SBarry Smith 96740244768SBarry Smith Options Database: 96840244768SBarry Smith + -snes_fas_levels - The number of levels 96940244768SBarry Smith . -snes_fas_cycles<1> - The number of cycles -- 1 for V, 2 for W 97040244768SBarry Smith . -snes_fas_type<additive,multiplicative,full,kaskade> - Additive or multiplicative cycle 97140244768SBarry Smith . -snes_fas_galerkin<PETSC_FALSE> - Form coarse problems by projection back upon the fine problem 97240244768SBarry Smith . -snes_fas_smoothup<1> - The number of iterations of the post-smoother 97340244768SBarry Smith . -snes_fas_smoothdown<1> - The number of iterations of the pre-smoother 97440244768SBarry Smith . -snes_fas_monitor - Monitor progress of all of the levels 97540244768SBarry Smith . -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS 97640244768SBarry Smith . -fas_levels_snes_ - SNES options for all smoothers 97740244768SBarry Smith . -fas_levels_cycle_snes_ - SNES options for all cycles 97840244768SBarry Smith . -fas_levels_i_snes_ - SNES options for the smoothers on level i 97940244768SBarry Smith . -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i 98040244768SBarry Smith - -fas_coarse_snes_ - SNES options for the coarsest smoother 98140244768SBarry Smith 98240244768SBarry Smith Notes: 98340244768SBarry Smith The organization of the FAS solver is slightly different from the organization of PCMG 98440244768SBarry Smith As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance. 98540244768SBarry Smith The cycle SNES instance may be used for monitoring convergence on a particular level. 98640244768SBarry Smith 98740244768SBarry Smith Level: beginner 98840244768SBarry Smith 98940244768SBarry Smith References: 99040244768SBarry Smith . 1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers", 99140244768SBarry Smith SIAM Review, 57(4), 2015 99240244768SBarry Smith 99340244768SBarry Smith .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types) 99440244768SBarry Smith M*/ 99540244768SBarry Smith 99640244768SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes) 99740244768SBarry Smith { 99840244768SBarry Smith SNES_FAS *fas; 99940244768SBarry Smith PetscErrorCode ierr; 100040244768SBarry Smith 100140244768SBarry Smith PetscFunctionBegin; 100240244768SBarry Smith snes->ops->destroy = SNESDestroy_FAS; 100340244768SBarry Smith snes->ops->setup = SNESSetUp_FAS; 100440244768SBarry Smith snes->ops->setfromoptions = SNESSetFromOptions_FAS; 100540244768SBarry Smith snes->ops->view = SNESView_FAS; 100640244768SBarry Smith snes->ops->solve = SNESSolve_FAS; 100740244768SBarry Smith snes->ops->reset = SNESReset_FAS; 100840244768SBarry Smith 100940244768SBarry Smith snes->usesksp = PETSC_FALSE; 1010efd4aadfSBarry Smith snes->usesnpc = PETSC_FALSE; 101140244768SBarry Smith 101240244768SBarry Smith if (!snes->tolerancesset) { 101340244768SBarry Smith snes->max_funcs = 30000; 101440244768SBarry Smith snes->max_its = 10000; 101540244768SBarry Smith } 101640244768SBarry Smith 10174fc747eaSLawrence Mitchell snes->alwayscomputesfinalresidual = PETSC_TRUE; 10184fc747eaSLawrence Mitchell 101940244768SBarry Smith ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr); 102040244768SBarry Smith 102140244768SBarry Smith snes->data = (void*) fas; 102240244768SBarry Smith fas->level = 0; 102340244768SBarry Smith fas->levels = 1; 102440244768SBarry Smith fas->n_cycles = 1; 102540244768SBarry Smith fas->max_up_it = 1; 102640244768SBarry Smith fas->max_down_it = 1; 102740244768SBarry Smith fas->smoothu = NULL; 102840244768SBarry Smith fas->smoothd = NULL; 102940244768SBarry Smith fas->next = NULL; 103040244768SBarry Smith fas->previous = NULL; 103140244768SBarry Smith fas->fine = snes; 103240244768SBarry Smith fas->interpolate = NULL; 103340244768SBarry Smith fas->restrct = NULL; 103440244768SBarry Smith fas->inject = NULL; 103540244768SBarry Smith fas->usedmfornumberoflevels = PETSC_FALSE; 103640244768SBarry Smith fas->fastype = SNES_FAS_MULTIPLICATIVE; 103740244768SBarry Smith fas->full_downsweep = PETSC_FALSE; 1038*6dfbafefSToby Isaac fas->full_total = PETSC_FALSE; 103940244768SBarry Smith 104040244768SBarry Smith fas->eventsmoothsetup = 0; 104140244768SBarry Smith fas->eventsmoothsolve = 0; 104240244768SBarry Smith fas->eventresidual = 0; 104340244768SBarry Smith fas->eventinterprestrict = 0; 104440244768SBarry Smith PetscFunctionReturn(0); 104540244768SBarry Smith } 1046