xref: /petsc/src/snes/impls/fas/fas.c (revision ba67282182bbff18a5c2e8153be2769e671cfaf7)
1421d9b32SPeter Brune /* Defines the basic SNES object */
2a055c8caSBarry Smith #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnes.h"  I*/
3421d9b32SPeter Brune 
434d65b3cSPeter Brune const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","SNESFASType","SNES_FAS",0};
507144faaSPeter Brune 
640244768SBarry Smith static PetscErrorCode SNESReset_FAS(SNES snes)
7421d9b32SPeter Brune {
877df8cc4SPeter Brune   PetscErrorCode ierr  = 0;
9421d9b32SPeter Brune   SNES_FAS       * fas = (SNES_FAS*)snes->data;
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   if (fas->galerkin) {
197cd33a7bSPeter Brune     ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr);
207cd33a7bSPeter Brune     ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr);
217cd33a7bSPeter Brune   }
22d477d801SPeter Brune   if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);}
23421d9b32SPeter Brune   PetscFunctionReturn(0);
24421d9b32SPeter Brune }
25421d9b32SPeter Brune 
2640244768SBarry Smith static PetscErrorCode SNESDestroy_FAS(SNES snes)
27421d9b32SPeter Brune {
28421d9b32SPeter Brune   SNES_FAS       * fas = (SNES_FAS*)snes->data;
29742fe5e2SPeter Brune   PetscErrorCode ierr  = 0;
30421d9b32SPeter Brune 
31421d9b32SPeter Brune   PetscFunctionBegin;
32421d9b32SPeter Brune   /* recursively resets and then destroys */
3379d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
341aa26658SKarl Rupp   if (fas->next) {
351aa26658SKarl Rupp     ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
361aa26658SKarl Rupp   }
37421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
38421d9b32SPeter Brune   PetscFunctionReturn(0);
39421d9b32SPeter Brune }
40421d9b32SPeter Brune 
4140244768SBarry Smith static PetscErrorCode SNESSetUp_FAS(SNES snes)
42421d9b32SPeter Brune {
4348bfdf8aSPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
44421d9b32SPeter Brune   PetscErrorCode ierr;
45d1adcc6fSPeter Brune   PetscInt       dm_levels;
463dccd265SPeter Brune   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
47ab8d36c9SPeter Brune   SNES           next;
480a7266b2SPatrick Farrell   PetscBool      isFine, hasCreateRestriction;
49f89ba88eSPeter Brune   SNESLineSearch linesearch;
50f89ba88eSPeter Brune   SNESLineSearch slinesearch;
51f89ba88eSPeter Brune   void           *lsprectx,*lspostctx;
526b2b7091SBarry Smith   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
536b2b7091SBarry Smith   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
54eff52c0eSPeter Brune 
556b2b7091SBarry Smith   PetscFunctionBegin;
56ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
57ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
58d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
59d1adcc6fSPeter Brune     dm_levels++;
60cc05f883SPeter Brune     if (dm_levels > fas->levels) {
612e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
623dccd265SPeter Brune       vec_sol              = snes->vec_sol;
633dccd265SPeter Brune       vec_func             = snes->vec_func;
643dccd265SPeter Brune       vec_sol_update       = snes->vec_sol_update;
653dccd265SPeter Brune       vec_rhs              = snes->vec_rhs;
660298fd71SBarry Smith       snes->vec_sol        = NULL;
670298fd71SBarry Smith       snes->vec_func       = NULL;
680298fd71SBarry Smith       snes->vec_sol_update = NULL;
690298fd71SBarry Smith       snes->vec_rhs        = NULL;
703dccd265SPeter Brune 
713dccd265SPeter Brune       /* reset the number of levels */
720298fd71SBarry Smith       ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr);
73cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
743dccd265SPeter Brune 
753dccd265SPeter Brune       snes->vec_sol        = vec_sol;
763dccd265SPeter Brune       snes->vec_func       = vec_func;
773dccd265SPeter Brune       snes->vec_rhs        = vec_rhs;
783dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
79d1adcc6fSPeter Brune     }
80d1adcc6fSPeter Brune   }
81ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
82ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
833dccd265SPeter Brune 
84fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
85cc05f883SPeter Brune 
86ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
87ab8d36c9SPeter Brune   if (!fas->smoothd) {
88ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
89ab8d36c9SPeter Brune   }
90ab8d36c9SPeter Brune 
9179d9a41aSPeter Brune   if (snes->dm) {
92ab8d36c9SPeter Brune     /* set the smoother DMs properly */
93ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
94ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
9579d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
96ab8d36c9SPeter Brune     if (next) {
9779d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
98ab8d36c9SPeter Brune       if (!next->dm) {
99ce94432eSBarry Smith         ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr);
100ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
10179d9a41aSPeter Brune       }
10279d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
10379d9a41aSPeter Brune       if (!fas->interpolate) {
104ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
105bccf9bb3SJed Brown         if (!fas->restrct) {
1060a7266b2SPatrick Farrell           ierr = DMHasCreateRestriction(next->dm, &hasCreateRestriction);CHKERRQ(ierr);
1070a7266b2SPatrick Farrell           /* DM can create restrictions, use that */
1080a7266b2SPatrick Farrell           if (hasCreateRestriction) {
1090a7266b2SPatrick Farrell             ierr = DMCreateRestriction(next->dm, snes->dm, &fas->restrct);CHKERRQ(ierr);
1100a7266b2SPatrick Farrell           } else {
111bccf9bb3SJed Brown             ierr         = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
11279d9a41aSPeter Brune             fas->restrct = fas->interpolate;
11379d9a41aSPeter Brune           }
114bccf9bb3SJed Brown         }
1150a7266b2SPatrick Farrell       }
11679d9a41aSPeter Brune       /* set the injection from the DM */
11779d9a41aSPeter Brune       if (!fas->inject) {
11823e68893SLawrence Mitchell         PetscBool flg;
11923e68893SLawrence Mitchell         ierr = DMHasCreateInjection(next->dm, &flg);CHKERRQ(ierr);
12023e68893SLawrence Mitchell         if (flg) {
1216dbf9973SLawrence Mitchell           ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr);
12279d9a41aSPeter Brune         }
12379d9a41aSPeter Brune       }
12479d9a41aSPeter Brune     }
12523e68893SLawrence Mitchell   }
12679d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
12779d9a41aSPeter Brune   if (fas->galerkin) {
1281aa26658SKarl Rupp     if (next) {
12925acbd8eSLisandro Dalcin       ierr = SNESSetFunction(next, NULL, SNESFASGalerkinFunctionDefault, next);CHKERRQ(ierr);
1301aa26658SKarl Rupp     }
1311aa26658SKarl Rupp     if (fas->smoothd && fas->level != fas->levels - 1) {
13225acbd8eSLisandro Dalcin       ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinFunctionDefault, snes);CHKERRQ(ierr);
1331aa26658SKarl Rupp     }
1341aa26658SKarl Rupp     if (fas->smoothu && fas->level != fas->levels - 1) {
13525acbd8eSLisandro Dalcin       ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinFunctionDefault, snes);CHKERRQ(ierr);
1361aa26658SKarl Rupp     }
13779d9a41aSPeter Brune   }
13879d9a41aSPeter Brune 
139534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
140534ebe21SPeter Brune   if (fas->smoothd) {
141bc3f2f05SPeter Brune     if (fas->level == 0 && fas->levels != 1) {
142365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
143534ebe21SPeter Brune     } else {
144365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
145534ebe21SPeter Brune     }
1467fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr);
147534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
1487601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
1497601faf0SJed Brown     ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr);
1506b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
1516b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
1526b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
1536b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
154f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
1550dd27c6cSPeter Brune 
1560dd27c6cSPeter Brune     fas->smoothd->vec_sol        = snes->vec_sol;
1570dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
1580dd27c6cSPeter Brune     fas->smoothd->vec_sol_update = snes->vec_sol_update;
1590dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
1600dd27c6cSPeter Brune     fas->smoothd->vec_func       = snes->vec_func;
1610dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
1620dd27c6cSPeter Brune 
163217044c2SLisandro Dalcin     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,fas->smoothd,0,0,0);CHKERRQ(ierr);}
1640dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr);
165217044c2SLisandro Dalcin     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,fas->smoothd,0,0,0);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     }
1757fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr);
176534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
1777601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
1787601faf0SJed Brown     ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr);
1796b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
1806b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
1816b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
1826b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
183f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
1840dd27c6cSPeter Brune 
1850dd27c6cSPeter Brune     fas->smoothu->vec_sol        = snes->vec_sol;
1860dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
1870dd27c6cSPeter Brune     fas->smoothu->vec_sol_update = snes->vec_sol_update;
1880dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
1890dd27c6cSPeter Brune     fas->smoothu->vec_func       = snes->vec_func;
1900dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
1910dd27c6cSPeter Brune 
192217044c2SLisandro Dalcin     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,fas->smoothu,0,0,0);CHKERRQ(ierr);}
1930dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr);
194217044c2SLisandro Dalcin     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,fas->smoothu,0,0,0);CHKERRQ(ierr);}
1950dd27c6cSPeter Brune 
196534ebe21SPeter Brune   }
197d06165b7SPeter Brune 
198ab8d36c9SPeter Brune   if (next) {
19979d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
200ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
201ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
2027fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr);
2037601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
2047601faf0SJed Brown     ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr);
2056b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2066b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2076b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2086b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
209f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
210ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
21179d9a41aSPeter Brune   }
2126273346dSPeter Brune   /* setup FAS work vectors */
2136273346dSPeter Brune   if (fas->galerkin) {
2146273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
2156273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
2166273346dSPeter Brune   }
217421d9b32SPeter Brune   PetscFunctionReturn(0);
218421d9b32SPeter Brune }
219421d9b32SPeter Brune 
22040244768SBarry Smith static PetscErrorCode SNESSetFromOptions_FAS(PetscOptionItems *PetscOptionsObject,SNES snes)
221421d9b32SPeter Brune {
222ee78dd50SPeter Brune   SNES_FAS       *fas   = (SNES_FAS*) snes->data;
223ee78dd50SPeter Brune   PetscInt       levels = 1;
22487f44e3fSPeter Brune   PetscBool      flg    = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE,continuationflg = PETSC_FALSE;
225421d9b32SPeter Brune   PetscErrorCode ierr;
22607144faaSPeter Brune   SNESFASType    fastype;
227fde0ff24SPeter Brune   const char     *optionsprefix;
228f1c6b773SPeter Brune   SNESLineSearch linesearch;
22966585501SPeter Brune   PetscInt       m, n_up, n_down;
230ab8d36c9SPeter Brune   SNES           next;
231ab8d36c9SPeter Brune   PetscBool      isFine;
232421d9b32SPeter Brune 
233421d9b32SPeter Brune   PetscFunctionBegin;
234ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
235e55864a3SBarry Smith   ierr = PetscOptionsHead(PetscOptionsObject,"SNESFAS Options-----------------------------------");CHKERRQ(ierr);
236ee78dd50SPeter Brune 
237ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
238ab8d36c9SPeter Brune   if (isFine) {
239ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
240c732cbdbSBarry Smith     if (!flg && snes->dm) {
241c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
242c732cbdbSBarry Smith       levels++;
243d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
244c732cbdbSBarry Smith     }
2450298fd71SBarry Smith     ierr    = SNESFASSetLevels(snes, levels, NULL);CHKERRQ(ierr);
24607144faaSPeter Brune     fastype = fas->fastype;
24707144faaSPeter Brune     ierr    = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
24807144faaSPeter Brune     if (flg) {
24907144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
25007144faaSPeter Brune     }
251ee78dd50SPeter Brune 
252fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
253ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
254ab8d36c9SPeter Brune     if (flg) {
255ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
256fde0ff24SPeter Brune     }
25787f44e3fSPeter Brune     ierr = PetscOptionsBool("-snes_fas_continuation","Corrected grid-sequence continuation","SNESFASSetContinuation",fas->continuation,&continuationflg,&flg);CHKERRQ(ierr);
25887f44e3fSPeter Brune     if (flg) {
25987f44e3fSPeter Brune       ierr = SNESFASSetContinuation(snes,continuationflg);CHKERRQ(ierr);
26087f44e3fSPeter Brune     }
261fde0ff24SPeter Brune 
262ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
263ab8d36c9SPeter Brune     if (flg) {
264ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
265ab8d36c9SPeter Brune     }
266ee78dd50SPeter Brune 
267928e959bSPeter Brune     if (fas->fastype == SNES_FAS_FULL) {
268*ba672821SBarry 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);
269*ba672821SBarry Smith       if (flg) {ierr = SNESFASFullSetDownSweep(snes,fas->full_downsweep);CHKERRQ(ierr);}
270928e959bSPeter Brune     }
271928e959bSPeter Brune 
27266585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
273162d76ddSPeter Brune 
27466585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
275162d76ddSPeter Brune 
276d142ab34SLawrence Mitchell     {
277d142ab34SLawrence Mitchell       PetscViewer viewer;
278d142ab34SLawrence Mitchell       PetscViewerFormat format;
27916413a6aSBarry Smith       ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject) snes)->options,((PetscObject)snes)->prefix,"-snes_fas_monitor",&viewer,&format,&monflg);CHKERRQ(ierr);
280d142ab34SLawrence Mitchell       if (monflg) {
281d142ab34SLawrence Mitchell         PetscViewerAndFormat *vf;
282d142ab34SLawrence Mitchell         ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
283d142ab34SLawrence Mitchell         ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
284d142ab34SLawrence Mitchell         ierr = SNESFASSetMonitor(snes,vf,PETSC_TRUE);CHKERRQ(ierr);
285d142ab34SLawrence Mitchell       }
286d142ab34SLawrence Mitchell     }
2870dd27c6cSPeter Brune     flg    = PETSC_FALSE;
2880dd27c6cSPeter Brune     monflg = PETSC_TRUE;
2890dd27c6cSPeter Brune     ierr   = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr);
2900dd27c6cSPeter Brune     if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);}
291ab8d36c9SPeter Brune   }
292ee78dd50SPeter Brune 
293421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
2948cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
295162d76ddSPeter Brune   if (upflg) {
29666585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
297162d76ddSPeter Brune   }
298162d76ddSPeter Brune   if (downflg) {
29966585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
300162d76ddSPeter Brune   }
301eff52c0eSPeter Brune 
3029e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
3039e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
3049e764e56SPeter Brune     if (!snes->linesearch) {
3057601faf0SJed Brown       ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
3061a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
3079e764e56SPeter Brune     }
3089e764e56SPeter Brune   }
3099e764e56SPeter Brune 
310ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
311ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
312ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
313421d9b32SPeter Brune   PetscFunctionReturn(0);
314421d9b32SPeter Brune }
315421d9b32SPeter Brune 
3169804daf3SBarry Smith #include <petscdraw.h>
31740244768SBarry Smith static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
318421d9b32SPeter Brune {
319421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
320656ede7eSPeter Brune   PetscBool      isFine,iascii,isdraw;
321ab8d36c9SPeter Brune   PetscInt       i;
322421d9b32SPeter Brune   PetscErrorCode ierr;
323ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
324421d9b32SPeter Brune 
325421d9b32SPeter Brune   PetscFunctionBegin;
326ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
327ab8d36c9SPeter Brune   if (isFine) {
328251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
329656ede7eSPeter Brune     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
330421d9b32SPeter Brune     if (iascii) {
331efd4aadfSBarry Smith       ierr = PetscViewerASCIIPrintf(viewer, "  type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
332ab8d36c9SPeter Brune       if (fas->galerkin) {
333ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
334421d9b32SPeter Brune       } else {
335ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"  Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
336421d9b32SPeter Brune       }
337ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
338ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
339ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
340ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
341ab8d36c9SPeter Brune         if (!i) {
342ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
343421d9b32SPeter Brune         } else {
344ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
345421d9b32SPeter Brune         }
346ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
347166b3ea4SJed Brown         if (smoothd) {
348ab8d36c9SPeter Brune           ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
349166b3ea4SJed Brown         } else {
350166b3ea4SJed Brown           ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
351166b3ea4SJed Brown         }
352ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
353ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
354ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
355ab8d36c9SPeter Brune         } else if (i) {
356ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"  Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
357ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
358166b3ea4SJed Brown           if (smoothu) {
359ab8d36c9SPeter Brune             ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
360166b3ea4SJed Brown           } else {
361166b3ea4SJed Brown             ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
362166b3ea4SJed Brown           }
363ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
364ab8d36c9SPeter Brune         }
365ab8d36c9SPeter Brune       }
366656ede7eSPeter Brune     } else if (isdraw) {
367656ede7eSPeter Brune       PetscDraw draw;
368b4375e8dSPeter Brune       PetscReal x,w,y,bottom,th,wth;
369656ede7eSPeter Brune       SNES_FAS  *curfas = fas;
370656ede7eSPeter Brune       ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
371656ede7eSPeter Brune       ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
372656ede7eSPeter Brune       ierr   = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr);
373656ede7eSPeter Brune       bottom = y - th;
374656ede7eSPeter Brune       while (curfas) {
375b4375e8dSPeter Brune         if (!curfas->smoothu) {
376656ede7eSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
377656ede7eSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
378656ede7eSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
379b4375e8dSPeter Brune         } else {
380b4375e8dSPeter Brune           w    = 0.5*PetscMin(1.0-x,x);
381b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
382b4375e8dSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
383b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
384b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
385b4375e8dSPeter Brune           if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr);
386b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
387b4375e8dSPeter Brune         }
388656ede7eSPeter Brune         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
389656ede7eSPeter Brune         bottom -= 5*th;
3901aa26658SKarl Rupp         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
3910298fd71SBarry Smith         else curfas = NULL;
392656ede7eSPeter Brune       }
393421d9b32SPeter Brune     }
394ab8d36c9SPeter Brune   }
395421d9b32SPeter Brune   PetscFunctionReturn(0);
396421d9b32SPeter Brune }
397421d9b32SPeter Brune 
39839bd7f45SPeter Brune /*
39939bd7f45SPeter Brune Defines the action of the downsmoother
40039bd7f45SPeter Brune  */
40140244768SBarry Smith static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
402b9c2fdf1SPeter Brune {
40339bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
404742fe5e2SPeter Brune   SNESConvergedReason reason;
405ab8d36c9SPeter Brune   Vec                 FPC;
406ab8d36c9SPeter Brune   SNES                smoothd;
4076cbb2f26SLawrence Mitchell   PetscBool           flg;
4080dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
4096e111a19SKarl Rupp 
410421d9b32SPeter Brune   PetscFunctionBegin;
411ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
412e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
413217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,smoothd,B,X,0);CHKERRQ(ierr);}
414ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
415217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,smoothd,B,X,0);CHKERRQ(ierr);}
416742fe5e2SPeter Brune   /* check convergence reason for the smoother */
417ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
4181ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
419742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
420742fe5e2SPeter Brune     PetscFunctionReturn(0);
421742fe5e2SPeter Brune   }
4226cbb2f26SLawrence Mitchell 
4230298fd71SBarry Smith   ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr);
4246cbb2f26SLawrence Mitchell   ierr = SNESGetAlwaysComputesFinalResidual(smoothd, &flg);CHKERRQ(ierr);
4256cbb2f26SLawrence Mitchell   if (!flg) {
4266cbb2f26SLawrence Mitchell     ierr = SNESComputeFunction(smoothd, X, FPC);CHKERRQ(ierr);
4276cbb2f26SLawrence Mitchell   }
4284b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
42971dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
43039bd7f45SPeter Brune   PetscFunctionReturn(0);
43139bd7f45SPeter Brune }
43239bd7f45SPeter Brune 
43339bd7f45SPeter Brune 
43439bd7f45SPeter Brune /*
43507144faaSPeter Brune Defines the action of the upsmoother
43639bd7f45SPeter Brune  */
43740244768SBarry Smith static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
4380adebc6cSBarry Smith {
43939bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
44039bd7f45SPeter Brune   SNESConvergedReason reason;
441ab8d36c9SPeter Brune   Vec                 FPC;
442ab8d36c9SPeter Brune   SNES                smoothu;
4436cbb2f26SLawrence Mitchell   PetscBool           flg;
4440dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
445ab8d36c9SPeter Brune 
4466e111a19SKarl Rupp   PetscFunctionBegin;
447ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
448217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,smoothu,0,0,0);CHKERRQ(ierr);}
449ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
450217044c2SLisandro Dalcin   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,smoothu,0,0,0);CHKERRQ(ierr);}
45139bd7f45SPeter Brune   /* check convergence reason for the smoother */
452ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
4531ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
45439bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
45539bd7f45SPeter Brune     PetscFunctionReturn(0);
45639bd7f45SPeter Brune   }
4570298fd71SBarry Smith   ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr);
4586cbb2f26SLawrence Mitchell   ierr = SNESGetAlwaysComputesFinalResidual(smoothu, &flg);CHKERRQ(ierr);
4596cbb2f26SLawrence Mitchell   if (!flg) {
4606cbb2f26SLawrence Mitchell     ierr = SNESComputeFunction(smoothu, X, FPC);CHKERRQ(ierr);
4616cbb2f26SLawrence Mitchell   }
4624b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
46371dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
46439bd7f45SPeter Brune   PetscFunctionReturn(0);
46539bd7f45SPeter Brune }
46639bd7f45SPeter Brune 
467938e4a01SJed Brown /*@
468938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
469938e4a01SJed Brown 
470938e4a01SJed Brown    Collective
471938e4a01SJed Brown 
472938e4a01SJed Brown    Input Arguments:
473938e4a01SJed Brown .  snes - SNESFAS
474938e4a01SJed Brown 
475938e4a01SJed Brown    Output Arguments:
476938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
477938e4a01SJed Brown 
478938e4a01SJed Brown    Level: developer
479938e4a01SJed Brown 
480938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
481938e4a01SJed Brown @*/
482938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
483938e4a01SJed Brown {
484938e4a01SJed Brown   PetscErrorCode ierr;
485938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
486938e4a01SJed Brown 
487938e4a01SJed Brown   PetscFunctionBegin;
4881aa26658SKarl Rupp   if (fas->rscale) {
4891aa26658SKarl Rupp     ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);
490f5af7f23SKarl Rupp   } else if (fas->inject) {
4912a7a6963SBarry Smith     ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr);
49213903a91SSatish Balay   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");
493938e4a01SJed Brown   PetscFunctionReturn(0);
494938e4a01SJed Brown }
495938e4a01SJed Brown 
496e9923e8dSJed Brown /*@
497e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
498e9923e8dSJed Brown 
499e9923e8dSJed Brown    Collective
500e9923e8dSJed Brown 
501e9923e8dSJed Brown    Input Arguments:
502e9923e8dSJed Brown +  fine - SNES from which to restrict
503e9923e8dSJed Brown -  Xfine - vector to restrict
504e9923e8dSJed Brown 
505e9923e8dSJed Brown    Output Arguments:
506e9923e8dSJed Brown .  Xcoarse - result of restriction
507e9923e8dSJed Brown 
508e9923e8dSJed Brown    Level: developer
509e9923e8dSJed Brown 
510e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
511e9923e8dSJed Brown @*/
512e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
513e9923e8dSJed Brown {
514e9923e8dSJed Brown   PetscErrorCode ierr;
515e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
516e9923e8dSJed Brown 
517e9923e8dSJed Brown   PetscFunctionBegin;
518e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
519e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
520e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
521e9923e8dSJed Brown   if (fas->inject) {
522e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
523e9923e8dSJed Brown   } else {
524e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
525e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
526e9923e8dSJed Brown   }
527e9923e8dSJed Brown   PetscFunctionReturn(0);
528e9923e8dSJed Brown }
529e9923e8dSJed Brown 
53039bd7f45SPeter Brune /*
53139bd7f45SPeter Brune 
53239bd7f45SPeter Brune Performs the FAS coarse correction as:
53339bd7f45SPeter Brune 
534b5527d98SMatthew G. Knepley fine problem:   F(x) = b
535b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c
53639bd7f45SPeter Brune 
537b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
53839bd7f45SPeter Brune 
53939bd7f45SPeter Brune  */
5400adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
5410adebc6cSBarry Smith {
54239bd7f45SPeter Brune   PetscErrorCode      ierr;
54339bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
54439bd7f45SPeter Brune   SNESConvergedReason reason;
545ab8d36c9SPeter Brune   SNES                next;
546ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
5470dd27c6cSPeter Brune   SNES_FAS            *fasc;
5485fd66863SKarl Rupp 
54939bd7f45SPeter Brune   PetscFunctionBegin;
550ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
551ab8d36c9SPeter Brune   if (next) {
5520dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
5530dd27c6cSPeter Brune 
554ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
555ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
556ab8d36c9SPeter Brune 
557ab8d36c9SPeter Brune     X_c  = next->vec_sol;
558ab8d36c9SPeter Brune     Xo_c = next->work[0];
559ab8d36c9SPeter Brune     F_c  = next->vec_func;
560ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
561efe1f98aSPeter Brune 
562217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
563938e4a01SJed Brown     ierr = SNESFASRestrict(snes, X, Xo_c);CHKERRQ(ierr);
5645609cbf2SMatthew G. Knepley     /* restrict the defect: R(F(x) - b) */
565ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
566217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
5670dd27c6cSPeter Brune 
568217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
5695609cbf2SMatthew G. Knepley     /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
570ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
571217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
5720dd27c6cSPeter Brune 
5730dd27c6cSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
574e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
575b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
576e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
577ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
578ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
579c90fad12SPeter Brune 
580c90fad12SPeter Brune     /* recurse to the next level */
581e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
582ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
583ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
584742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
585742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
586742fe5e2SPeter Brune       PetscFunctionReturn(0);
587742fe5e2SPeter Brune     }
588fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
589fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
5900dd27c6cSPeter Brune 
591217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
592ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
593fb8dc43dSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) X_c, "Coarse correction");CHKERRQ(ierr);
594fb8dc43dSMatthew G. Knepley     ierr = VecViewFromOptions(X_c, NULL, "-fas_coarse_solution_view");CHKERRQ(ierr);
595fb8dc43dSMatthew G. Knepley     ierr = PetscObjectSetName((PetscObject) X_new, "Updated Fine solution");CHKERRQ(ierr);
596fb8dc43dSMatthew G. Knepley     ierr = VecViewFromOptions(X_new, NULL, "-fas_levels_1_solution_view");CHKERRQ(ierr);
597217044c2SLisandro Dalcin     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,snes,0,0,0);CHKERRQ(ierr);}
598293a7e31SPeter Brune   }
59939bd7f45SPeter Brune   PetscFunctionReturn(0);
60039bd7f45SPeter Brune }
60139bd7f45SPeter Brune 
60239bd7f45SPeter Brune /*
60339bd7f45SPeter Brune 
60439bd7f45SPeter Brune The additive cycle looks like:
60539bd7f45SPeter Brune 
60607144faaSPeter Brune xhat = x
60707144faaSPeter Brune xhat = dS(x, b)
60807144faaSPeter Brune x = coarsecorrection(xhat, b_d)
60907144faaSPeter Brune x = x + nu*(xhat - x);
61039bd7f45SPeter Brune (optional) x = uS(x, b)
61139bd7f45SPeter Brune 
61239bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
61339bd7f45SPeter Brune 
61439bd7f45SPeter Brune  */
61540244768SBarry Smith static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
6160adebc6cSBarry Smith {
61707144faaSPeter Brune   Vec                  F, B, Xhat;
61822c1e704SPeter Brune   Vec                  X_c, Xo_c, F_c, B_c;
61939bd7f45SPeter Brune   PetscErrorCode       ierr;
62007144faaSPeter Brune   SNESConvergedReason  reason;
62122c1e704SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
622422a814eSBarry Smith   SNESLineSearchReason lsresult;
623ab8d36c9SPeter Brune   SNES                 next;
624ab8d36c9SPeter Brune   Mat                  restrct, interpolate;
6250dd27c6cSPeter Brune   SNES_FAS             *fas = (SNES_FAS*)snes->data,*fasc;
6260dd27c6cSPeter Brune 
62739bd7f45SPeter Brune   PetscFunctionBegin;
628ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
62939bd7f45SPeter Brune   F    = snes->vec_func;
63039bd7f45SPeter Brune   B    = snes->vec_rhs;
631e7f468e7SPeter Brune   Xhat = snes->work[1];
63207144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
63307144faaSPeter Brune   /* recurse first */
634ab8d36c9SPeter Brune   if (next) {
6350dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
636ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
637ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
638217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
63907144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
640217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
641c2a02606SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
642ab8d36c9SPeter Brune     X_c  = next->vec_sol;
643ab8d36c9SPeter Brune     Xo_c = next->work[0];
644ab8d36c9SPeter Brune     F_c  = next->vec_func;
645ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
64639bd7f45SPeter Brune 
647938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
64807144faaSPeter Brune     /* restrict the defect */
649ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
65007144faaSPeter Brune 
65107144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
652217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
653ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
654217044c2SLisandro Dalcin     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,next,0,0,0);CHKERRQ(ierr);}
655e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
656b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
657e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
65807144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
65907144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
66007144faaSPeter Brune 
66107144faaSPeter Brune     /* recurse */
662e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
663ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
66407144faaSPeter Brune 
66507144faaSPeter Brune     /* smooth on this level */
66691f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
66707144faaSPeter Brune 
668ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
66907144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
67007144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
67107144faaSPeter Brune       PetscFunctionReturn(0);
67207144faaSPeter Brune     }
67307144faaSPeter Brune 
67407144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
675c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
676ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
67707144faaSPeter Brune 
678ddebd997SPeter Brune     /* additive correction of the coarse direction*/
679f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
680422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr);
681422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
682422a814eSBarry Smith     if (lsresult) {
6839e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6849e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6859e764e56SPeter Brune         PetscFunctionReturn(0);
6869e764e56SPeter Brune       }
6879e764e56SPeter Brune     }
68807144faaSPeter Brune   } else {
68991f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
69007144faaSPeter Brune   }
69139bd7f45SPeter Brune   PetscFunctionReturn(0);
69239bd7f45SPeter Brune }
69339bd7f45SPeter Brune 
69439bd7f45SPeter Brune /*
69539bd7f45SPeter Brune 
69639bd7f45SPeter Brune Defines the FAS cycle as:
69739bd7f45SPeter Brune 
698b5527d98SMatthew G. Knepley fine problem: F(x) = b
69939bd7f45SPeter Brune coarse problem: F^c(x) = b^c
70039bd7f45SPeter Brune 
701b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
70239bd7f45SPeter Brune 
70339bd7f45SPeter Brune correction:
70439bd7f45SPeter Brune 
70539bd7f45SPeter Brune x = x + I(x^c - Rx)
70639bd7f45SPeter Brune 
70739bd7f45SPeter Brune  */
70840244768SBarry Smith static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
7090adebc6cSBarry Smith {
71039bd7f45SPeter Brune 
71139bd7f45SPeter Brune   PetscErrorCode ierr;
71239bd7f45SPeter Brune   Vec            F,B;
71334d65b3cSPeter Brune   SNES           next;
71439bd7f45SPeter Brune 
71539bd7f45SPeter Brune   PetscFunctionBegin;
71639bd7f45SPeter Brune   F = snes->vec_func;
71739bd7f45SPeter Brune   B = snes->vec_rhs;
71839bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
71934d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
72091f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
72134d65b3cSPeter Brune   if (next) {
7228c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
72391f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
724fe6f9142SPeter Brune   }
725fa9694d7SPeter Brune   PetscFunctionReturn(0);
726421d9b32SPeter Brune }
727421d9b32SPeter Brune 
72840244768SBarry Smith static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
7298c94862eSPeter Brune {
7308c94862eSPeter Brune   SNES           next;
7318c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7328c94862eSPeter Brune   PetscBool      isFine;
7338c94862eSPeter Brune   PetscErrorCode ierr;
7348c94862eSPeter Brune 
7358c94862eSPeter Brune   PetscFunctionBegin;
7368c94862eSPeter Brune   /* pre-smooth -- just update using the pre-smoother */
7378c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
7388c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
7398c94862eSPeter Brune   fas->full_stage = 0;
7408c94862eSPeter Brune   if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);}
7418c94862eSPeter Brune   PetscFunctionReturn(0);
7428c94862eSPeter Brune }
7438c94862eSPeter Brune 
74440244768SBarry Smith static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
7458c94862eSPeter Brune {
7468c94862eSPeter Brune   PetscErrorCode ierr;
7478c94862eSPeter Brune   Vec            F,B;
7488c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7498c94862eSPeter Brune   PetscBool      isFine;
7508c94862eSPeter Brune   SNES           next;
7518c94862eSPeter Brune 
7528c94862eSPeter Brune   PetscFunctionBegin;
7538c94862eSPeter Brune   F = snes->vec_func;
7548c94862eSPeter Brune   B = snes->vec_rhs;
7558c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
7568c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
7578c94862eSPeter Brune 
7588c94862eSPeter Brune   if (isFine) {
7598c94862eSPeter Brune     ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr);
7608c94862eSPeter Brune   }
7618c94862eSPeter Brune 
7628c94862eSPeter Brune   if (fas->full_stage == 0) {
763928e959bSPeter Brune     /* downsweep */
7648c94862eSPeter Brune     if (next) {
7658c94862eSPeter Brune       if (fas->level != 1) next->max_its += 1;
766e140b65aSPatrick Farrell       if (fas->full_downsweep) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
767e140b65aSPatrick Farrell       fas->full_downsweep = PETSC_TRUE;
7688c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
7698c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7708c94862eSPeter Brune       if (fas->level != 1) next->max_its -= 1;
7718c94862eSPeter Brune     } else {
772a3a80b83SMatthew G. Knepley       /* The smoother on the coarse level is the coarse solver */
7738c94862eSPeter Brune       ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7748c94862eSPeter Brune     }
7758c94862eSPeter Brune     fas->full_stage = 1;
7768c94862eSPeter Brune   } else if (fas->full_stage == 1) {
7778c94862eSPeter Brune     if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
7788c94862eSPeter Brune     if (next) {
7798c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
7808c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7818c94862eSPeter Brune     }
7828c94862eSPeter Brune   }
7838c94862eSPeter Brune   /* final v-cycle */
7848c94862eSPeter Brune   if (isFine) {
7858c94862eSPeter Brune     if (next) {
7868c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
7878c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7888c94862eSPeter Brune     }
7898c94862eSPeter Brune   }
7908c94862eSPeter Brune   PetscFunctionReturn(0);
7918c94862eSPeter Brune }
7928c94862eSPeter Brune 
79340244768SBarry Smith static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
79434d65b3cSPeter Brune {
79534d65b3cSPeter Brune   PetscErrorCode ierr;
79634d65b3cSPeter Brune   Vec            F,B;
79734d65b3cSPeter Brune   SNES           next;
79834d65b3cSPeter Brune 
79934d65b3cSPeter Brune   PetscFunctionBegin;
80034d65b3cSPeter Brune   F = snes->vec_func;
80134d65b3cSPeter Brune   B = snes->vec_rhs;
80234d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
80334d65b3cSPeter Brune   if (next) {
80434d65b3cSPeter Brune     ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
80534d65b3cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
80634d65b3cSPeter Brune   } else {
80734d65b3cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
80834d65b3cSPeter Brune   }
80934d65b3cSPeter Brune   PetscFunctionReturn(0);
81034d65b3cSPeter Brune }
81134d65b3cSPeter Brune 
812fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE;
813fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n"
814fffbeea8SBarry Smith                             "  title = {Composing Scalable Nonlinear Algebraic Solvers},\n"
815fffbeea8SBarry Smith                             "  author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n"
816fffbeea8SBarry Smith                             "  year = 2013,\n"
817fffbeea8SBarry Smith                             "  type = Preprint,\n"
818fffbeea8SBarry Smith                             "  number = {ANL/MCS-P2010-0112},\n"
819fffbeea8SBarry Smith                             "  institution = {Argonne National Laboratory}\n}\n";
820fffbeea8SBarry Smith 
82140244768SBarry Smith static PetscErrorCode SNESSolve_FAS(SNES snes)
822421d9b32SPeter Brune {
823fa9694d7SPeter Brune   PetscErrorCode ierr;
824fe6f9142SPeter Brune   PetscInt       i, maxits;
825ddb5aff1SPeter Brune   Vec            X, F;
826fe6f9142SPeter Brune   PetscReal      fnorm;
827b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
828b17ce1afSJed Brown   DM             dm;
829e70c42e5SPeter Brune   PetscBool      isFine;
830b17ce1afSJed Brown 
831421d9b32SPeter Brune   PetscFunctionBegin;
832c579b300SPatrick Farrell 
8336c4ed002SBarry 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);
834c579b300SPatrick Farrell 
835fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
836fe6f9142SPeter Brune   maxits       = snes->max_its;      /* maximum number of iterations */
837fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
838fa9694d7SPeter Brune   X            = snes->vec_sol;
839f5a6d4f9SBarry Smith   F            = snes->vec_func;
840293a7e31SPeter Brune 
84118a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
842293a7e31SPeter Brune   /*norm setup */
843e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
844fe6f9142SPeter Brune   snes->iter = 0;
845fe6f9142SPeter Brune   snes->norm = 0.;
846e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
847e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
848217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
849fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
850217044c2SLisandro Dalcin     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,snes,0,0,0);CHKERRQ(ierr);}
8511aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
852e4ed7901SPeter Brune 
853fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
854422a814eSBarry Smith   SNESCheckFunctionNorm(snes,fnorm);
855e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
856fe6f9142SPeter Brune   snes->norm = fnorm;
857e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
858a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
859fe6f9142SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
860fe6f9142SPeter Brune 
861fe6f9142SPeter Brune   /* test convergence */
862fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
863fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
864e4ed7901SPeter Brune 
865b17ce1afSJed Brown 
866b9c2fdf1SPeter Brune   if (isFine) {
867b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
868b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
869b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
870b17ce1afSJed Brown       DM dmcoarse;
871b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
872b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
873b17ce1afSJed Brown       dm   = dmcoarse;
874b17ce1afSJed Brown     }
875b9c2fdf1SPeter Brune   }
876b17ce1afSJed Brown 
877fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
878fe6f9142SPeter Brune     /* Call general purpose update function */
879646217ecSPeter Brune 
880fe6f9142SPeter Brune     if (snes->ops->update) {
881fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
882fe6f9142SPeter Brune     }
88307144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
88491f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
8858c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
88691f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
8878c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_FULL) {
8888c94862eSPeter Brune       ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr);
88934d65b3cSPeter Brune     } else if (fas->fastype ==SNES_FAS_KASKADE) {
89034d65b3cSPeter Brune       ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr);
8916c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");
892742fe5e2SPeter Brune 
893742fe5e2SPeter Brune     /* check for FAS cycle divergence */
8941aa26658SKarl Rupp     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0);
895b9c2fdf1SPeter Brune 
896c90fad12SPeter Brune     /* Monitor convergence */
897e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
898c90fad12SPeter Brune     snes->iter = i+1;
899e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
900a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
901c90fad12SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
902c90fad12SPeter Brune     /* Test for convergence */
90366585501SPeter Brune     if (isFine) {
904b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
905c90fad12SPeter Brune       if (snes->reason) break;
906fe6f9142SPeter Brune     }
90766585501SPeter Brune   }
908fe6f9142SPeter Brune   if (i == maxits) {
909fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
910fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
911fe6f9142SPeter Brune   }
912421d9b32SPeter Brune   PetscFunctionReturn(0);
913421d9b32SPeter Brune }
91440244768SBarry Smith 
91540244768SBarry Smith /*MC
91640244768SBarry Smith 
91740244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
91840244768SBarry Smith 
91940244768SBarry Smith    The nonlinear problem is solved by correction using coarse versions
92040244768SBarry Smith    of the nonlinear problem.  This problem is perturbed so that a projected
92140244768SBarry Smith    solution of the fine problem elicits no correction from the coarse problem.
92240244768SBarry Smith 
92340244768SBarry Smith Options Database:
92440244768SBarry Smith +   -snes_fas_levels -  The number of levels
92540244768SBarry Smith .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
92640244768SBarry Smith .   -snes_fas_type<additive,multiplicative,full,kaskade>  -  Additive or multiplicative cycle
92740244768SBarry Smith .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
92840244768SBarry Smith .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
92940244768SBarry Smith .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
93040244768SBarry Smith .   -snes_fas_monitor -  Monitor progress of all of the levels
93140244768SBarry Smith .   -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS
93240244768SBarry Smith .   -fas_levels_snes_ -  SNES options for all smoothers
93340244768SBarry Smith .   -fas_levels_cycle_snes_ -  SNES options for all cycles
93440244768SBarry Smith .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
93540244768SBarry Smith .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
93640244768SBarry Smith -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
93740244768SBarry Smith 
93840244768SBarry Smith Notes:
93940244768SBarry Smith    The organization of the FAS solver is slightly different from the organization of PCMG
94040244768SBarry Smith    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
94140244768SBarry Smith    The cycle SNES instance may be used for monitoring convergence on a particular level.
94240244768SBarry Smith 
94340244768SBarry Smith Level: beginner
94440244768SBarry Smith 
94540244768SBarry Smith    References:
94640244768SBarry Smith . 1. -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
94740244768SBarry Smith    SIAM Review, 57(4), 2015
94840244768SBarry Smith 
94940244768SBarry Smith .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
95040244768SBarry Smith M*/
95140244768SBarry Smith 
95240244768SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
95340244768SBarry Smith {
95440244768SBarry Smith   SNES_FAS       *fas;
95540244768SBarry Smith   PetscErrorCode ierr;
95640244768SBarry Smith 
95740244768SBarry Smith   PetscFunctionBegin;
95840244768SBarry Smith   snes->ops->destroy        = SNESDestroy_FAS;
95940244768SBarry Smith   snes->ops->setup          = SNESSetUp_FAS;
96040244768SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
96140244768SBarry Smith   snes->ops->view           = SNESView_FAS;
96240244768SBarry Smith   snes->ops->solve          = SNESSolve_FAS;
96340244768SBarry Smith   snes->ops->reset          = SNESReset_FAS;
96440244768SBarry Smith 
96540244768SBarry Smith   snes->usesksp = PETSC_FALSE;
966efd4aadfSBarry Smith   snes->usesnpc = PETSC_FALSE;
96740244768SBarry Smith 
96840244768SBarry Smith   if (!snes->tolerancesset) {
96940244768SBarry Smith     snes->max_funcs = 30000;
97040244768SBarry Smith     snes->max_its   = 10000;
97140244768SBarry Smith   }
97240244768SBarry Smith 
9734fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
9744fc747eaSLawrence Mitchell 
97540244768SBarry Smith   ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr);
97640244768SBarry Smith 
97740244768SBarry Smith   snes->data                  = (void*) fas;
97840244768SBarry Smith   fas->level                  = 0;
97940244768SBarry Smith   fas->levels                 = 1;
98040244768SBarry Smith   fas->n_cycles               = 1;
98140244768SBarry Smith   fas->max_up_it              = 1;
98240244768SBarry Smith   fas->max_down_it            = 1;
98340244768SBarry Smith   fas->smoothu                = NULL;
98440244768SBarry Smith   fas->smoothd                = NULL;
98540244768SBarry Smith   fas->next                   = NULL;
98640244768SBarry Smith   fas->previous               = NULL;
98740244768SBarry Smith   fas->fine                   = snes;
98840244768SBarry Smith   fas->interpolate            = NULL;
98940244768SBarry Smith   fas->restrct                = NULL;
99040244768SBarry Smith   fas->inject                 = NULL;
99140244768SBarry Smith   fas->usedmfornumberoflevels = PETSC_FALSE;
99240244768SBarry Smith   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
99340244768SBarry Smith   fas->full_downsweep         = PETSC_FALSE;
99440244768SBarry Smith 
99540244768SBarry Smith   fas->eventsmoothsetup    = 0;
99640244768SBarry Smith   fas->eventsmoothsolve    = 0;
99740244768SBarry Smith   fas->eventresidual       = 0;
99840244768SBarry Smith   fas->eventinterprestrict = 0;
99940244768SBarry Smith   PetscFunctionReturn(0);
100040244768SBarry Smith }
1001