xref: /petsc/src/snes/impls/fas/fas.c (revision 6cbb2f269e2ef07e969434a4c788decb1633cd2c)
1421d9b32SPeter Brune /* Defines the basic SNES object */
2a055c8caSBarry Smith #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnes.h"  I*/
3421d9b32SPeter Brune 
434d65b3cSPeter Brune const char *const SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","FULL","KASKADE","SNESFASType","SNES_FAS",0};
507144faaSPeter Brune 
66273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES,Vec,Vec,void*);
7efe1f98aSPeter Brune 
8421d9b32SPeter Brune #undef __FUNCT__
9421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
1040244768SBarry Smith static PetscErrorCode SNESReset_FAS(SNES snes)
11421d9b32SPeter Brune {
1277df8cc4SPeter Brune   PetscErrorCode ierr  = 0;
13421d9b32SPeter Brune   SNES_FAS       * fas = (SNES_FAS*)snes->data;
14421d9b32SPeter Brune 
15421d9b32SPeter Brune   PetscFunctionBegin;
16ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
17ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
183dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
19bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
20bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
21bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
227cd33a7bSPeter Brune   if (fas->galerkin) {
237cd33a7bSPeter Brune     ierr = VecDestroy(&fas->Xg);CHKERRQ(ierr);
247cd33a7bSPeter Brune     ierr = VecDestroy(&fas->Fg);CHKERRQ(ierr);
257cd33a7bSPeter Brune   }
26d477d801SPeter Brune   if (fas->next) {ierr = SNESReset(fas->next);CHKERRQ(ierr);}
27421d9b32SPeter Brune   PetscFunctionReturn(0);
28421d9b32SPeter Brune }
29421d9b32SPeter Brune 
30421d9b32SPeter Brune #undef __FUNCT__
31421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
3240244768SBarry Smith static PetscErrorCode SNESDestroy_FAS(SNES snes)
33421d9b32SPeter Brune {
34421d9b32SPeter Brune   SNES_FAS       * fas = (SNES_FAS*)snes->data;
35742fe5e2SPeter Brune   PetscErrorCode ierr  = 0;
36421d9b32SPeter Brune 
37421d9b32SPeter Brune   PetscFunctionBegin;
38421d9b32SPeter Brune   /* recursively resets and then destroys */
3979d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
401aa26658SKarl Rupp   if (fas->next) {
411aa26658SKarl Rupp     ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
421aa26658SKarl Rupp   }
43421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
44421d9b32SPeter Brune   PetscFunctionReturn(0);
45421d9b32SPeter Brune }
46421d9b32SPeter Brune 
47421d9b32SPeter Brune #undef __FUNCT__
48421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
4940244768SBarry Smith static PetscErrorCode SNESSetUp_FAS(SNES snes)
50421d9b32SPeter Brune {
5148bfdf8aSPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
52421d9b32SPeter Brune   PetscErrorCode ierr;
53d1adcc6fSPeter Brune   PetscInt       dm_levels;
543dccd265SPeter Brune   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
55ab8d36c9SPeter Brune   SNES           next;
56ab8d36c9SPeter Brune   PetscBool      isFine;
57f89ba88eSPeter Brune   SNESLineSearch linesearch;
58f89ba88eSPeter Brune   SNESLineSearch slinesearch;
59f89ba88eSPeter Brune   void           *lsprectx,*lspostctx;
606b2b7091SBarry Smith   PetscErrorCode (*precheck)(SNESLineSearch,Vec,Vec,PetscBool*,void*);
616b2b7091SBarry Smith   PetscErrorCode (*postcheck)(SNESLineSearch,Vec,Vec,Vec,PetscBool*,PetscBool*,void*);
62eff52c0eSPeter Brune 
636b2b7091SBarry Smith   PetscFunctionBegin;
64ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
65ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
66d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
67d1adcc6fSPeter Brune     dm_levels++;
68cc05f883SPeter Brune     if (dm_levels > fas->levels) {
692e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
703dccd265SPeter Brune       vec_sol              = snes->vec_sol;
713dccd265SPeter Brune       vec_func             = snes->vec_func;
723dccd265SPeter Brune       vec_sol_update       = snes->vec_sol_update;
733dccd265SPeter Brune       vec_rhs              = snes->vec_rhs;
740298fd71SBarry Smith       snes->vec_sol        = NULL;
750298fd71SBarry Smith       snes->vec_func       = NULL;
760298fd71SBarry Smith       snes->vec_sol_update = NULL;
770298fd71SBarry Smith       snes->vec_rhs        = NULL;
783dccd265SPeter Brune 
793dccd265SPeter Brune       /* reset the number of levels */
800298fd71SBarry Smith       ierr = SNESFASSetLevels(snes,dm_levels,NULL);CHKERRQ(ierr);
81cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
823dccd265SPeter Brune 
833dccd265SPeter Brune       snes->vec_sol        = vec_sol;
843dccd265SPeter Brune       snes->vec_func       = vec_func;
853dccd265SPeter Brune       snes->vec_rhs        = vec_rhs;
863dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
87d1adcc6fSPeter Brune     }
88d1adcc6fSPeter Brune   }
89ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
90ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
913dccd265SPeter Brune 
92fa0ddf94SBarry Smith   ierr = SNESSetWorkVecs(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
93cc05f883SPeter Brune 
94ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
95ab8d36c9SPeter Brune   if (!fas->smoothd) {
96ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
97ab8d36c9SPeter Brune   }
98ab8d36c9SPeter Brune 
9979d9a41aSPeter Brune   if (snes->dm) {
100ab8d36c9SPeter Brune     /* set the smoother DMs properly */
101ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
102ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
10379d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
104ab8d36c9SPeter Brune     if (next) {
10579d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
106ab8d36c9SPeter Brune       if (!next->dm) {
107ce94432eSBarry Smith         ierr = DMCoarsen(snes->dm, PetscObjectComm((PetscObject)next), &next->dm);CHKERRQ(ierr);
108ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
10979d9a41aSPeter Brune       }
11079d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
11179d9a41aSPeter Brune       if (!fas->interpolate) {
112ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
113bccf9bb3SJed Brown         if (!fas->restrct) {
114bccf9bb3SJed Brown           ierr         = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
11579d9a41aSPeter Brune           fas->restrct = fas->interpolate;
11679d9a41aSPeter Brune         }
117bccf9bb3SJed Brown       }
11879d9a41aSPeter Brune       /* set the injection from the DM */
11979d9a41aSPeter Brune       if (!fas->inject) {
1206dbf9973SLawrence Mitchell         ierr = DMCreateInjection(next->dm, snes->dm, &fas->inject);CHKERRQ(ierr);
12179d9a41aSPeter Brune       }
12279d9a41aSPeter Brune     }
12379d9a41aSPeter Brune   }
12479d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
12579d9a41aSPeter Brune   if (fas->galerkin) {
1261aa26658SKarl Rupp     if (next) {
1270298fd71SBarry Smith       ierr = SNESSetFunction(next, NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
1281aa26658SKarl Rupp     }
1291aa26658SKarl Rupp     if (fas->smoothd && fas->level != fas->levels - 1) {
1300298fd71SBarry Smith       ierr = SNESSetFunction(fas->smoothd, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
1311aa26658SKarl Rupp     }
1321aa26658SKarl Rupp     if (fas->smoothu && fas->level != fas->levels - 1) {
1330298fd71SBarry Smith       ierr = SNESSetFunction(fas->smoothu, NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
1341aa26658SKarl Rupp     }
13579d9a41aSPeter Brune   }
13679d9a41aSPeter Brune 
137534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
138534ebe21SPeter Brune   if (fas->smoothd) {
139bc3f2f05SPeter Brune     if (fas->level == 0 && fas->levels != 1) {
140365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
141534ebe21SPeter Brune     } else {
142365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
143534ebe21SPeter Brune     }
1447fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothd);CHKERRQ(ierr);
145534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
1467601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
1477601faf0SJed Brown     ierr = SNESGetLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr);
1486b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
1496b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
1506b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
1516b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
152f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
1530dd27c6cSPeter Brune 
1540dd27c6cSPeter Brune     fas->smoothd->vec_sol        = snes->vec_sol;
1550dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
1560dd27c6cSPeter Brune     fas->smoothd->vec_sol_update = snes->vec_sol_update;
1570dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
1580dd27c6cSPeter Brune     fas->smoothd->vec_func       = snes->vec_func;
1590dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
1600dd27c6cSPeter Brune 
1610dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1620dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothd);CHKERRQ(ierr);
1630dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
164534ebe21SPeter Brune   }
165534ebe21SPeter Brune 
166534ebe21SPeter Brune   /* sets the up (post) smoother's default norm and sets it from options */
167534ebe21SPeter Brune   if (fas->smoothu) {
168534ebe21SPeter Brune     if (fas->level != fas->levels - 1) {
169365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr);
170534ebe21SPeter Brune     } else {
171365a6726SPeter Brune       ierr = SNESSetNormSchedule(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
172534ebe21SPeter Brune     }
1737fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)fas->smoothu);CHKERRQ(ierr);
174534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
1757601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
1767601faf0SJed Brown     ierr = SNESGetLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr);
1776b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
1786b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
1796b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
1806b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
181f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
1820dd27c6cSPeter Brune 
1830dd27c6cSPeter Brune     fas->smoothu->vec_sol        = snes->vec_sol;
1840dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol);CHKERRQ(ierr);
1850dd27c6cSPeter Brune     fas->smoothu->vec_sol_update = snes->vec_sol_update;
1860dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_sol_update);CHKERRQ(ierr);
1870dd27c6cSPeter Brune     fas->smoothu->vec_func       = snes->vec_func;
1880dd27c6cSPeter Brune     ierr                         = PetscObjectReference((PetscObject)snes->vec_func);CHKERRQ(ierr);
1890dd27c6cSPeter Brune 
1900dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventBegin(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1910dd27c6cSPeter Brune     ierr = SNESSetUp(fas->smoothu);CHKERRQ(ierr);
1920dd27c6cSPeter Brune     if (fas->eventsmoothsetup) {ierr = PetscLogEventEnd(fas->eventsmoothsetup,0,0,0,0);CHKERRQ(ierr);}
1930dd27c6cSPeter Brune 
194534ebe21SPeter Brune   }
195d06165b7SPeter Brune 
196ab8d36c9SPeter Brune   if (next) {
19779d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
198ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
199ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
2007fce8c19SPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)next);CHKERRQ(ierr);
2017601faf0SJed Brown     ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
2027601faf0SJed Brown     ierr = SNESGetLineSearch(fas->next,&slinesearch);CHKERRQ(ierr);
2036b2b7091SBarry Smith     ierr = SNESLineSearchGetPreCheck(linesearch,&precheck,&lsprectx);CHKERRQ(ierr);
2046b2b7091SBarry Smith     ierr = SNESLineSearchGetPostCheck(linesearch,&postcheck,&lspostctx);CHKERRQ(ierr);
2056b2b7091SBarry Smith     ierr = SNESLineSearchSetPreCheck(slinesearch,precheck,lsprectx);CHKERRQ(ierr);
2066b2b7091SBarry Smith     ierr = SNESLineSearchSetPostCheck(slinesearch,postcheck,lspostctx);CHKERRQ(ierr);
207f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
208ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
20979d9a41aSPeter Brune   }
2106273346dSPeter Brune   /* setup FAS work vectors */
2116273346dSPeter Brune   if (fas->galerkin) {
2126273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
2136273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
2146273346dSPeter Brune   }
215421d9b32SPeter Brune   PetscFunctionReturn(0);
216421d9b32SPeter Brune }
217421d9b32SPeter Brune 
218421d9b32SPeter Brune #undef __FUNCT__
219421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
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) {
268928e959bSPeter Brune       ierr   = PetscOptionsBool("-snes_fas_full_downsweep","Smooth on the initial upsweep for full FAS cycles","SNESFASFullSetDownSweep",fas->full_downsweep,&fas->full_downsweep,&flg);CHKERRQ(ierr);
269928e959bSPeter Brune       if (flg) {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;
279d142ab34SLawrence Mitchell       ierr = PetscOptionsGetViewer(PetscObjectComm((PetscObject)snes),((PetscObject)snes)->prefix,
280d142ab34SLawrence Mitchell                                    "-snes_fas_monitor",&viewer,&format,&monflg);CHKERRQ(ierr);
281d142ab34SLawrence Mitchell       if (monflg) {
282d142ab34SLawrence Mitchell         PetscViewerAndFormat *vf;
283d142ab34SLawrence Mitchell         ierr = PetscViewerAndFormatCreate(viewer,format,&vf);CHKERRQ(ierr);
284d142ab34SLawrence Mitchell         ierr = PetscObjectDereference((PetscObject)viewer);CHKERRQ(ierr);
285d142ab34SLawrence Mitchell         ierr = SNESFASSetMonitor(snes,vf,PETSC_TRUE);CHKERRQ(ierr);
286d142ab34SLawrence Mitchell       }
287d142ab34SLawrence Mitchell     }
2880dd27c6cSPeter Brune     flg    = PETSC_FALSE;
2890dd27c6cSPeter Brune     monflg = PETSC_TRUE;
2900dd27c6cSPeter Brune     ierr   = PetscOptionsBool("-snes_fas_log","Log times for each FAS level","SNESFASSetLog",monflg,&monflg,&flg);CHKERRQ(ierr);
2910dd27c6cSPeter Brune     if (flg) {ierr = SNESFASSetLog(snes,monflg);CHKERRQ(ierr);}
292ab8d36c9SPeter Brune   }
293ee78dd50SPeter Brune 
294421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
2958cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
296162d76ddSPeter Brune   if (upflg) {
29766585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
298162d76ddSPeter Brune   }
299162d76ddSPeter Brune   if (downflg) {
30066585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
301162d76ddSPeter Brune   }
302eff52c0eSPeter Brune 
3039e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
3049e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
3059e764e56SPeter Brune     if (!snes->linesearch) {
3067601faf0SJed Brown       ierr = SNESGetLineSearch(snes, &linesearch);CHKERRQ(ierr);
3071a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
3089e764e56SPeter Brune     }
3099e764e56SPeter Brune   }
3109e764e56SPeter Brune 
311ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
312ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
313ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
314421d9b32SPeter Brune   PetscFunctionReturn(0);
315421d9b32SPeter Brune }
316421d9b32SPeter Brune 
3179804daf3SBarry Smith #include <petscdraw.h>
318421d9b32SPeter Brune #undef __FUNCT__
319421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
32040244768SBarry Smith static PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
321421d9b32SPeter Brune {
322421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS*) snes->data;
323656ede7eSPeter Brune   PetscBool      isFine,iascii,isdraw;
324ab8d36c9SPeter Brune   PetscInt       i;
325421d9b32SPeter Brune   PetscErrorCode ierr;
326ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
327421d9b32SPeter Brune 
328421d9b32SPeter Brune   PetscFunctionBegin;
329ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
330ab8d36c9SPeter Brune   if (isFine) {
331251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr);
332656ede7eSPeter Brune     ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERDRAW,&isdraw);CHKERRQ(ierr);
333421d9b32SPeter Brune     if (iascii) {
334ab8d36c9SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
335ab8d36c9SPeter Brune       if (fas->galerkin) {
336ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
337421d9b32SPeter Brune       } else {
338ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
339421d9b32SPeter Brune       }
340ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
341ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
342ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
343ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
344ab8d36c9SPeter Brune         if (!i) {
345ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
346421d9b32SPeter Brune         } else {
347ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
348421d9b32SPeter Brune         }
349ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
350166b3ea4SJed Brown         if (smoothd) {
351ab8d36c9SPeter Brune           ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
352166b3ea4SJed Brown         } else {
353166b3ea4SJed Brown           ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
354166b3ea4SJed Brown         }
355ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
356ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
357ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
358ab8d36c9SPeter Brune         } else if (i) {
359ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
360ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
361166b3ea4SJed Brown           if (smoothu) {
362ab8d36c9SPeter Brune             ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
363166b3ea4SJed Brown           } else {
364166b3ea4SJed Brown             ierr = PetscViewerASCIIPrintf(viewer,"Not yet available\n");CHKERRQ(ierr);
365166b3ea4SJed Brown           }
366ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
367ab8d36c9SPeter Brune         }
368ab8d36c9SPeter Brune       }
369656ede7eSPeter Brune     } else if (isdraw) {
370656ede7eSPeter Brune       PetscDraw draw;
371b4375e8dSPeter Brune       PetscReal x,w,y,bottom,th,wth;
372656ede7eSPeter Brune       SNES_FAS  *curfas = fas;
373656ede7eSPeter Brune       ierr   = PetscViewerDrawGetDraw(viewer,0,&draw);CHKERRQ(ierr);
374656ede7eSPeter Brune       ierr   = PetscDrawGetCurrentPoint(draw,&x,&y);CHKERRQ(ierr);
375656ede7eSPeter Brune       ierr   = PetscDrawStringGetSize(draw,&wth,&th);CHKERRQ(ierr);
376656ede7eSPeter Brune       bottom = y - th;
377656ede7eSPeter Brune       while (curfas) {
378b4375e8dSPeter Brune         if (!curfas->smoothu) {
379656ede7eSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x,bottom);CHKERRQ(ierr);
380656ede7eSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
381656ede7eSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
382b4375e8dSPeter Brune         } else {
383b4375e8dSPeter Brune           w    = 0.5*PetscMin(1.0-x,x);
384b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x-w,bottom);CHKERRQ(ierr);
385b4375e8dSPeter Brune           if (curfas->smoothd) ierr = SNESView(curfas->smoothd,viewer);CHKERRQ(ierr);
386b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
387b4375e8dSPeter Brune           ierr = PetscDrawPushCurrentPoint(draw,x+w,bottom);CHKERRQ(ierr);
388b4375e8dSPeter Brune           if (curfas->smoothu) ierr = SNESView(curfas->smoothu,viewer);CHKERRQ(ierr);
389b4375e8dSPeter Brune           ierr = PetscDrawPopCurrentPoint(draw);CHKERRQ(ierr);
390b4375e8dSPeter Brune         }
391656ede7eSPeter Brune         /* this is totally bogus but we have no way of knowing how low the previous one was draw to */
392656ede7eSPeter Brune         bottom -= 5*th;
3931aa26658SKarl Rupp         if (curfas->next) curfas = (SNES_FAS*)curfas->next->data;
3940298fd71SBarry Smith         else curfas = NULL;
395656ede7eSPeter Brune       }
396421d9b32SPeter Brune     }
397ab8d36c9SPeter Brune   }
398421d9b32SPeter Brune   PetscFunctionReturn(0);
399421d9b32SPeter Brune }
400421d9b32SPeter Brune 
401421d9b32SPeter Brune #undef __FUNCT__
40291f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private"
40339bd7f45SPeter Brune /*
40439bd7f45SPeter Brune Defines the action of the downsmoother
40539bd7f45SPeter Brune  */
40640244768SBarry Smith static PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
407b9c2fdf1SPeter Brune {
40839bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
409742fe5e2SPeter Brune   SNESConvergedReason reason;
410ab8d36c9SPeter Brune   Vec                 FPC;
411ab8d36c9SPeter Brune   SNES                smoothd;
412*6cbb2f26SLawrence Mitchell   PetscBool           flg;
4130dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
4146e111a19SKarl Rupp 
415421d9b32SPeter Brune   PetscFunctionBegin;
416ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
417e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
4180dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
419ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
4200dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
421742fe5e2SPeter Brune   /* check convergence reason for the smoother */
422ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
4231ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
424742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
425742fe5e2SPeter Brune     PetscFunctionReturn(0);
426742fe5e2SPeter Brune   }
427*6cbb2f26SLawrence Mitchell 
4280298fd71SBarry Smith   ierr = SNESGetFunction(smoothd, &FPC, NULL, NULL);CHKERRQ(ierr);
429*6cbb2f26SLawrence Mitchell   ierr = SNESGetAlwaysComputesFinalResidual(smoothd, &flg);CHKERRQ(ierr);
430*6cbb2f26SLawrence Mitchell   if (!flg) {
431*6cbb2f26SLawrence Mitchell     ierr = SNESComputeFunction(smoothd, X, FPC);CHKERRQ(ierr);
432*6cbb2f26SLawrence Mitchell   }
4334b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
43471dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
43539bd7f45SPeter Brune   PetscFunctionReturn(0);
43639bd7f45SPeter Brune }
43739bd7f45SPeter Brune 
43839bd7f45SPeter Brune 
43939bd7f45SPeter Brune #undef __FUNCT__
44091f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private"
44139bd7f45SPeter Brune /*
44207144faaSPeter Brune Defines the action of the upsmoother
44339bd7f45SPeter Brune  */
44440244768SBarry Smith static PetscErrorCode SNESFASUpSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
4450adebc6cSBarry Smith {
44639bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
44739bd7f45SPeter Brune   SNESConvergedReason reason;
448ab8d36c9SPeter Brune   Vec                 FPC;
449ab8d36c9SPeter Brune   SNES                smoothu;
450*6cbb2f26SLawrence Mitchell   PetscBool           flg;
4510dd27c6cSPeter Brune   SNES_FAS            *fas = (SNES_FAS*) snes->data;
452ab8d36c9SPeter Brune 
4536e111a19SKarl Rupp   PetscFunctionBegin;
454ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
4550dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventBegin(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
456ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
4570dd27c6cSPeter Brune   if (fas->eventsmoothsolve) {ierr = PetscLogEventEnd(fas->eventsmoothsolve,0,0,0,0);CHKERRQ(ierr);}
45839bd7f45SPeter Brune   /* check convergence reason for the smoother */
459ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
4601ff805c3SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN || reason == SNES_DIVERGED_LINE_SEARCH)) {
46139bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
46239bd7f45SPeter Brune     PetscFunctionReturn(0);
46339bd7f45SPeter Brune   }
4640298fd71SBarry Smith   ierr = SNESGetFunction(smoothu, &FPC, NULL, NULL);CHKERRQ(ierr);
465*6cbb2f26SLawrence Mitchell   ierr = SNESGetAlwaysComputesFinalResidual(smoothu, &flg);CHKERRQ(ierr);
466*6cbb2f26SLawrence Mitchell   if (!flg) {
467*6cbb2f26SLawrence Mitchell     ierr = SNESComputeFunction(smoothu, X, FPC);CHKERRQ(ierr);
468*6cbb2f26SLawrence Mitchell   }
4694b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
47071dbe336SPeter Brune   if (fnorm) {ierr = VecNorm(F,NORM_2,fnorm);CHKERRQ(ierr);}
47139bd7f45SPeter Brune   PetscFunctionReturn(0);
47239bd7f45SPeter Brune }
47339bd7f45SPeter Brune 
47439bd7f45SPeter Brune #undef __FUNCT__
475938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
476938e4a01SJed Brown /*@
477938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
478938e4a01SJed Brown 
479938e4a01SJed Brown    Collective
480938e4a01SJed Brown 
481938e4a01SJed Brown    Input Arguments:
482938e4a01SJed Brown .  snes - SNESFAS
483938e4a01SJed Brown 
484938e4a01SJed Brown    Output Arguments:
485938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
486938e4a01SJed Brown 
487938e4a01SJed Brown    Level: developer
488938e4a01SJed Brown 
489938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
490938e4a01SJed Brown @*/
491938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
492938e4a01SJed Brown {
493938e4a01SJed Brown   PetscErrorCode ierr;
494938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
495938e4a01SJed Brown 
496938e4a01SJed Brown   PetscFunctionBegin;
4971aa26658SKarl Rupp   if (fas->rscale) {
4981aa26658SKarl Rupp     ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);
499f5af7f23SKarl Rupp   } else if (fas->inject) {
5002a7a6963SBarry Smith     ierr = MatCreateVecs(fas->inject,Xcoarse,NULL);CHKERRQ(ierr);
501ce94432eSBarry Smith   } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
502938e4a01SJed Brown   PetscFunctionReturn(0);
503938e4a01SJed Brown }
504938e4a01SJed Brown 
505e9923e8dSJed Brown #undef __FUNCT__
506e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
507e9923e8dSJed Brown /*@
508e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
509e9923e8dSJed Brown 
510e9923e8dSJed Brown    Collective
511e9923e8dSJed Brown 
512e9923e8dSJed Brown    Input Arguments:
513e9923e8dSJed Brown +  fine - SNES from which to restrict
514e9923e8dSJed Brown -  Xfine - vector to restrict
515e9923e8dSJed Brown 
516e9923e8dSJed Brown    Output Arguments:
517e9923e8dSJed Brown .  Xcoarse - result of restriction
518e9923e8dSJed Brown 
519e9923e8dSJed Brown    Level: developer
520e9923e8dSJed Brown 
521e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
522e9923e8dSJed Brown @*/
523e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
524e9923e8dSJed Brown {
525e9923e8dSJed Brown   PetscErrorCode ierr;
526e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
527e9923e8dSJed Brown 
528e9923e8dSJed Brown   PetscFunctionBegin;
529e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
530e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
531e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
532e9923e8dSJed Brown   if (fas->inject) {
533e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
534e9923e8dSJed Brown   } else {
535e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
536e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
537e9923e8dSJed Brown   }
538e9923e8dSJed Brown   PetscFunctionReturn(0);
539e9923e8dSJed Brown }
540e9923e8dSJed Brown 
541e9923e8dSJed Brown #undef __FUNCT__
5428c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection"
54339bd7f45SPeter Brune /*
54439bd7f45SPeter Brune 
54539bd7f45SPeter Brune Performs the FAS coarse correction as:
54639bd7f45SPeter Brune 
547b5527d98SMatthew G. Knepley fine problem:   F(x) = b
548b5527d98SMatthew G. Knepley coarse problem: F^c(x^c) = b^c
54939bd7f45SPeter Brune 
550b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
55139bd7f45SPeter Brune 
55239bd7f45SPeter Brune  */
5530adebc6cSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new)
5540adebc6cSBarry Smith {
55539bd7f45SPeter Brune   PetscErrorCode      ierr;
55639bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
55739bd7f45SPeter Brune   SNESConvergedReason reason;
558ab8d36c9SPeter Brune   SNES                next;
559ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
5600dd27c6cSPeter Brune   SNES_FAS            *fasc;
5615fd66863SKarl Rupp 
56239bd7f45SPeter Brune   PetscFunctionBegin;
563ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
564ab8d36c9SPeter Brune   if (next) {
5650dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
5660dd27c6cSPeter Brune 
567ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
568ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
569ab8d36c9SPeter Brune 
570ab8d36c9SPeter Brune     X_c  = next->vec_sol;
571ab8d36c9SPeter Brune     Xo_c = next->work[0];
572ab8d36c9SPeter Brune     F_c  = next->vec_func;
573ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
574efe1f98aSPeter Brune 
5750dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
576938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
5775609cbf2SMatthew G. Knepley     /* restrict the defect: R(F(x) - b) */
578ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
5790dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
5800dd27c6cSPeter Brune 
5810dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
5825609cbf2SMatthew G. Knepley     /* F_c = F^c(Rx) - R(F(x) - b) since the second term was sitting in next->vec_rhs */
583ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
5840dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
5850dd27c6cSPeter Brune 
5860dd27c6cSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
587e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
588b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
589e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
590ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
591ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
592c90fad12SPeter Brune 
593c90fad12SPeter Brune     /* recurse to the next level */
594e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
595ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
596ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
597742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
598742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
599742fe5e2SPeter Brune       PetscFunctionReturn(0);
600742fe5e2SPeter Brune     }
601fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
602fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
6030dd27c6cSPeter Brune 
6040dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventBegin(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
605ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
6060dd27c6cSPeter Brune     if (fasc->eventinterprestrict) {ierr = PetscLogEventEnd(fasc->eventinterprestrict,0,0,0,0);CHKERRQ(ierr);}
607293a7e31SPeter Brune   }
60839bd7f45SPeter Brune   PetscFunctionReturn(0);
60939bd7f45SPeter Brune }
61039bd7f45SPeter Brune 
61139bd7f45SPeter Brune #undef __FUNCT__
6122cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive"
61339bd7f45SPeter Brune /*
61439bd7f45SPeter Brune 
61539bd7f45SPeter Brune The additive cycle looks like:
61639bd7f45SPeter Brune 
61707144faaSPeter Brune xhat = x
61807144faaSPeter Brune xhat = dS(x, b)
61907144faaSPeter Brune x = coarsecorrection(xhat, b_d)
62007144faaSPeter Brune x = x + nu*(xhat - x);
62139bd7f45SPeter Brune (optional) x = uS(x, b)
62239bd7f45SPeter Brune 
62339bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
62439bd7f45SPeter Brune 
62539bd7f45SPeter Brune  */
62640244768SBarry Smith static PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X)
6270adebc6cSBarry Smith {
62807144faaSPeter Brune   Vec                  F, B, Xhat;
62922c1e704SPeter Brune   Vec                  X_c, Xo_c, F_c, B_c;
63039bd7f45SPeter Brune   PetscErrorCode       ierr;
63107144faaSPeter Brune   SNESConvergedReason  reason;
63222c1e704SPeter Brune   PetscReal            xnorm, fnorm, ynorm;
633422a814eSBarry Smith   SNESLineSearchReason lsresult;
634ab8d36c9SPeter Brune   SNES                 next;
635ab8d36c9SPeter Brune   Mat                  restrct, interpolate;
6360dd27c6cSPeter Brune   SNES_FAS             *fas = (SNES_FAS*)snes->data,*fasc;
6370dd27c6cSPeter Brune 
63839bd7f45SPeter Brune   PetscFunctionBegin;
639ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
64039bd7f45SPeter Brune   F    = snes->vec_func;
64139bd7f45SPeter Brune   B    = snes->vec_rhs;
642e7f468e7SPeter Brune   Xhat = snes->work[1];
64307144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
64407144faaSPeter Brune   /* recurse first */
645ab8d36c9SPeter Brune   if (next) {
6460dd27c6cSPeter Brune     fasc = (SNES_FAS*)next->data;
647ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
648ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
6490dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
65007144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
6510dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
652c2a02606SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
653ab8d36c9SPeter Brune     X_c  = next->vec_sol;
654ab8d36c9SPeter Brune     Xo_c = next->work[0];
655ab8d36c9SPeter Brune     F_c  = next->vec_func;
656ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
65739bd7f45SPeter Brune 
658938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
65907144faaSPeter Brune     /* restrict the defect */
660ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
66107144faaSPeter Brune 
66207144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
6630dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventBegin(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
664ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
6650dd27c6cSPeter Brune     if (fasc->eventresidual) {ierr = PetscLogEventEnd(fasc->eventresidual,0,0,0,0);CHKERRQ(ierr);}
666e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
667b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
668e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
66907144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
67007144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
67107144faaSPeter Brune 
67207144faaSPeter Brune     /* recurse */
673e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
674ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
67507144faaSPeter Brune 
67607144faaSPeter Brune     /* smooth on this level */
67791f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
67807144faaSPeter Brune 
679ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
68007144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
68107144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
68207144faaSPeter Brune       PetscFunctionReturn(0);
68307144faaSPeter Brune     }
68407144faaSPeter Brune 
68507144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
686c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
687ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
68807144faaSPeter Brune 
689ddebd997SPeter Brune     /* additive correction of the coarse direction*/
690f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
691422a814eSBarry Smith     ierr = SNESLineSearchGetReason(snes->linesearch, &lsresult);CHKERRQ(ierr);
692422a814eSBarry Smith     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
693422a814eSBarry Smith     if (lsresult) {
6949e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6959e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6969e764e56SPeter Brune         PetscFunctionReturn(0);
6979e764e56SPeter Brune       }
6989e764e56SPeter Brune     }
69907144faaSPeter Brune   } else {
70091f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
70107144faaSPeter Brune   }
70239bd7f45SPeter Brune   PetscFunctionReturn(0);
70339bd7f45SPeter Brune }
70439bd7f45SPeter Brune 
70539bd7f45SPeter Brune #undef __FUNCT__
7062cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative"
70739bd7f45SPeter Brune /*
70839bd7f45SPeter Brune 
70939bd7f45SPeter Brune Defines the FAS cycle as:
71039bd7f45SPeter Brune 
711b5527d98SMatthew G. Knepley fine problem: F(x) = b
71239bd7f45SPeter Brune coarse problem: F^c(x) = b^c
71339bd7f45SPeter Brune 
714b5527d98SMatthew G. Knepley b^c = F^c(Rx) - R(F(x) - b)
71539bd7f45SPeter Brune 
71639bd7f45SPeter Brune correction:
71739bd7f45SPeter Brune 
71839bd7f45SPeter Brune x = x + I(x^c - Rx)
71939bd7f45SPeter Brune 
72039bd7f45SPeter Brune  */
72140244768SBarry Smith static PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X)
7220adebc6cSBarry Smith {
72339bd7f45SPeter Brune 
72439bd7f45SPeter Brune   PetscErrorCode ierr;
72539bd7f45SPeter Brune   Vec            F,B;
72634d65b3cSPeter Brune   SNES           next;
72739bd7f45SPeter Brune 
72839bd7f45SPeter Brune   PetscFunctionBegin;
72939bd7f45SPeter Brune   F = snes->vec_func;
73039bd7f45SPeter Brune   B = snes->vec_rhs;
73139bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
73234d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
73391f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
73434d65b3cSPeter Brune   if (next) {
7358c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
73691f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
737fe6f9142SPeter Brune   }
738fa9694d7SPeter Brune   PetscFunctionReturn(0);
739421d9b32SPeter Brune }
740421d9b32SPeter Brune 
741421d9b32SPeter Brune #undef __FUNCT__
7428c94862eSPeter Brune #define __FUNCT__ "SNESFASCycleSetupPhase_Full"
74340244768SBarry Smith static PetscErrorCode SNESFASCycleSetupPhase_Full(SNES snes)
7448c94862eSPeter Brune {
7458c94862eSPeter Brune   SNES           next;
7468c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7478c94862eSPeter Brune   PetscBool      isFine;
7488c94862eSPeter Brune   PetscErrorCode ierr;
7498c94862eSPeter Brune 
7508c94862eSPeter Brune   PetscFunctionBegin;
7518c94862eSPeter Brune   /* pre-smooth -- just update using the pre-smoother */
7528c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
7538c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
7548c94862eSPeter Brune   fas->full_stage = 0;
7558c94862eSPeter Brune   if (next) {ierr = SNESFASCycleSetupPhase_Full(next);CHKERRQ(ierr);}
7568c94862eSPeter Brune   PetscFunctionReturn(0);
7578c94862eSPeter Brune }
7588c94862eSPeter Brune 
7598c94862eSPeter Brune #undef __FUNCT__
7608c94862eSPeter Brune #define __FUNCT__ "SNESFASCycle_Full"
76140244768SBarry Smith static PetscErrorCode SNESFASCycle_Full(SNES snes, Vec X)
7628c94862eSPeter Brune {
7638c94862eSPeter Brune   PetscErrorCode ierr;
7648c94862eSPeter Brune   Vec            F,B;
7658c94862eSPeter Brune   SNES_FAS       *fas = (SNES_FAS*)snes->data;
7668c94862eSPeter Brune   PetscBool      isFine;
7678c94862eSPeter Brune   SNES           next;
7688c94862eSPeter Brune 
7698c94862eSPeter Brune   PetscFunctionBegin;
7708c94862eSPeter Brune   F = snes->vec_func;
7718c94862eSPeter Brune   B = snes->vec_rhs;
7728c94862eSPeter Brune   ierr = SNESFASCycleIsFine(snes,&isFine);CHKERRQ(ierr);
7738c94862eSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
7748c94862eSPeter Brune 
7758c94862eSPeter Brune   if (isFine) {
7768c94862eSPeter Brune     ierr = SNESFASCycleSetupPhase_Full(snes);CHKERRQ(ierr);
7778c94862eSPeter Brune   }
7788c94862eSPeter Brune 
7798c94862eSPeter Brune   if (fas->full_stage == 0) {
780928e959bSPeter Brune     /* downsweep */
7818c94862eSPeter Brune     if (next) {
7828c94862eSPeter Brune       if (fas->level != 1) next->max_its += 1;
783928e959bSPeter Brune       if (fas->full_downsweep||isFine) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
7848c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
7858c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7868c94862eSPeter Brune       if (fas->level != 1) next->max_its -= 1;
7878c94862eSPeter Brune     } else {
788a3a80b83SMatthew G. Knepley       /* The smoother on the coarse level is the coarse solver */
7898c94862eSPeter Brune       ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7908c94862eSPeter Brune     }
7918c94862eSPeter Brune     fas->full_stage = 1;
7928c94862eSPeter Brune   } else if (fas->full_stage == 1) {
7938c94862eSPeter Brune     if (snes->iter == 0) {ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);}
7948c94862eSPeter Brune     if (next) {
7958c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
7968c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
7978c94862eSPeter Brune     }
7988c94862eSPeter Brune   }
7998c94862eSPeter Brune   /* final v-cycle */
8008c94862eSPeter Brune   if (isFine) {
8018c94862eSPeter Brune     if (next) {
8028c94862eSPeter Brune       ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
8038c94862eSPeter Brune       ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
8048c94862eSPeter Brune     }
8058c94862eSPeter Brune   }
8068c94862eSPeter Brune   PetscFunctionReturn(0);
8078c94862eSPeter Brune }
8088c94862eSPeter Brune 
8098c94862eSPeter Brune #undef __FUNCT__
81034d65b3cSPeter Brune #define __FUNCT__ "SNESFASCycle_Kaskade"
81140244768SBarry Smith static PetscErrorCode SNESFASCycle_Kaskade(SNES snes, Vec X)
81234d65b3cSPeter Brune {
81334d65b3cSPeter Brune   PetscErrorCode ierr;
81434d65b3cSPeter Brune   Vec            F,B;
81534d65b3cSPeter Brune   SNES           next;
81634d65b3cSPeter Brune 
81734d65b3cSPeter Brune   PetscFunctionBegin;
81834d65b3cSPeter Brune   F = snes->vec_func;
81934d65b3cSPeter Brune   B = snes->vec_rhs;
82034d65b3cSPeter Brune   ierr = SNESFASCycleGetCorrection(snes,&next);CHKERRQ(ierr);
82134d65b3cSPeter Brune   if (next) {
82234d65b3cSPeter Brune     ierr = SNESFASCoarseCorrection(snes,X,F,X);CHKERRQ(ierr);
82334d65b3cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
82434d65b3cSPeter Brune   } else {
82534d65b3cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes,B,X,F,&snes->norm);CHKERRQ(ierr);
82634d65b3cSPeter Brune   }
82734d65b3cSPeter Brune   PetscFunctionReturn(0);
82834d65b3cSPeter Brune }
82934d65b3cSPeter Brune 
830fffbeea8SBarry Smith PetscBool SNEScite = PETSC_FALSE;
831fffbeea8SBarry Smith const char SNESCitation[] = "@techreport{pbmkbsxt2012,\n"
832fffbeea8SBarry Smith                             "  title = {Composing Scalable Nonlinear Algebraic Solvers},\n"
833fffbeea8SBarry Smith                             "  author = {Peter Brune and Mathew Knepley and Barry Smith and Xuemin Tu},\n"
834fffbeea8SBarry Smith                             "  year = 2013,\n"
835fffbeea8SBarry Smith                             "  type = Preprint,\n"
836fffbeea8SBarry Smith                             "  number = {ANL/MCS-P2010-0112},\n"
837fffbeea8SBarry Smith                             "  institution = {Argonne National Laboratory}\n}\n";
838fffbeea8SBarry Smith 
83934d65b3cSPeter Brune #undef __FUNCT__
840421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
84140244768SBarry Smith static PetscErrorCode SNESSolve_FAS(SNES snes)
842421d9b32SPeter Brune {
843fa9694d7SPeter Brune   PetscErrorCode ierr;
844fe6f9142SPeter Brune   PetscInt       i, maxits;
845ddb5aff1SPeter Brune   Vec            X, F;
846fe6f9142SPeter Brune   PetscReal      fnorm;
847b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data,*ffas;
848b17ce1afSJed Brown   DM             dm;
849e70c42e5SPeter Brune   PetscBool      isFine;
850b17ce1afSJed Brown 
851421d9b32SPeter Brune   PetscFunctionBegin;
852c579b300SPatrick Farrell 
8536c4ed002SBarry 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);
854c579b300SPatrick Farrell 
855fffbeea8SBarry Smith   ierr = PetscCitationsRegister(SNESCitation,&SNEScite);CHKERRQ(ierr);
856fe6f9142SPeter Brune   maxits       = snes->max_its;      /* maximum number of iterations */
857fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
858fa9694d7SPeter Brune   X            = snes->vec_sol;
859f5a6d4f9SBarry Smith   F            = snes->vec_func;
860293a7e31SPeter Brune 
86118a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
862293a7e31SPeter Brune   /*norm setup */
863e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
864fe6f9142SPeter Brune   snes->iter = 0;
865fe6f9142SPeter Brune   snes->norm = 0.;
866e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
867e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
8680dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventBegin(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
869fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
8700dd27c6cSPeter Brune     if (fas->eventresidual) {ierr = PetscLogEventEnd(fas->eventresidual,0,0,0,0);CHKERRQ(ierr);}
8711aa26658SKarl Rupp   } else snes->vec_func_init_set = PETSC_FALSE;
872e4ed7901SPeter Brune 
873fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
874422a814eSBarry Smith   SNESCheckFunctionNorm(snes,fnorm);
875e04113cfSBarry Smith   ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
876fe6f9142SPeter Brune   snes->norm = fnorm;
877e04113cfSBarry Smith   ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
878a71f0d7dSBarry Smith   ierr       = SNESLogConvergenceHistory(snes,fnorm,0);CHKERRQ(ierr);
879fe6f9142SPeter Brune   ierr       = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
880fe6f9142SPeter Brune 
881fe6f9142SPeter Brune   /* test convergence */
882fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
883fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
884e4ed7901SPeter Brune 
885b17ce1afSJed Brown 
886b9c2fdf1SPeter Brune   if (isFine) {
887b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
888b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
889b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
890b17ce1afSJed Brown       DM dmcoarse;
891b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
892b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
893b17ce1afSJed Brown       dm   = dmcoarse;
894b17ce1afSJed Brown     }
895b9c2fdf1SPeter Brune   }
896b17ce1afSJed Brown 
897fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
898fe6f9142SPeter Brune     /* Call general purpose update function */
899646217ecSPeter Brune 
900fe6f9142SPeter Brune     if (snes->ops->update) {
901fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
902fe6f9142SPeter Brune     }
90307144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
90491f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
9058c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_ADDITIVE) {
90691f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
9078c94862eSPeter Brune     } else if (fas->fastype == SNES_FAS_FULL) {
9088c94862eSPeter Brune       ierr = SNESFASCycle_Full(snes, X);CHKERRQ(ierr);
90934d65b3cSPeter Brune     } else if (fas->fastype ==SNES_FAS_KASKADE) {
91034d65b3cSPeter Brune       ierr = SNESFASCycle_Kaskade(snes, X);CHKERRQ(ierr);
9116c4ed002SBarry Smith     } else SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE,"Unsupported FAS type");
912742fe5e2SPeter Brune 
913742fe5e2SPeter Brune     /* check for FAS cycle divergence */
9141aa26658SKarl Rupp     if (snes->reason != SNES_CONVERGED_ITERATING) PetscFunctionReturn(0);
915b9c2fdf1SPeter Brune 
916c90fad12SPeter Brune     /* Monitor convergence */
917e04113cfSBarry Smith     ierr       = PetscObjectSAWsTakeAccess((PetscObject)snes);CHKERRQ(ierr);
918c90fad12SPeter Brune     snes->iter = i+1;
919e04113cfSBarry Smith     ierr       = PetscObjectSAWsGrantAccess((PetscObject)snes);CHKERRQ(ierr);
920a71f0d7dSBarry Smith     ierr       = SNESLogConvergenceHistory(snes,snes->norm,0);CHKERRQ(ierr);
921c90fad12SPeter Brune     ierr       = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
922c90fad12SPeter Brune     /* Test for convergence */
92366585501SPeter Brune     if (isFine) {
924b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
925c90fad12SPeter Brune       if (snes->reason) break;
926fe6f9142SPeter Brune     }
92766585501SPeter Brune   }
928fe6f9142SPeter Brune   if (i == maxits) {
929fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
930fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
931fe6f9142SPeter Brune   }
932421d9b32SPeter Brune   PetscFunctionReturn(0);
933421d9b32SPeter Brune }
93440244768SBarry Smith 
93540244768SBarry Smith /*MC
93640244768SBarry Smith 
93740244768SBarry Smith SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
93840244768SBarry Smith 
93940244768SBarry Smith    The nonlinear problem is solved by correction using coarse versions
94040244768SBarry Smith    of the nonlinear problem.  This problem is perturbed so that a projected
94140244768SBarry Smith    solution of the fine problem elicits no correction from the coarse problem.
94240244768SBarry Smith 
94340244768SBarry Smith Options Database:
94440244768SBarry Smith +   -snes_fas_levels -  The number of levels
94540244768SBarry Smith .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
94640244768SBarry Smith .   -snes_fas_type<additive,multiplicative,full,kaskade>  -  Additive or multiplicative cycle
94740244768SBarry Smith .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
94840244768SBarry Smith .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
94940244768SBarry Smith .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
95040244768SBarry Smith .   -snes_fas_monitor -  Monitor progress of all of the levels
95140244768SBarry Smith .   -snes_fas_full_downsweep<PETSC_FALSE> - call the downsmooth on the initial downsweep of full FAS
95240244768SBarry Smith .   -fas_levels_snes_ -  SNES options for all smoothers
95340244768SBarry Smith .   -fas_levels_cycle_snes_ -  SNES options for all cycles
95440244768SBarry Smith .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
95540244768SBarry Smith .   -fas_levels_i_cycle_snes_ - SNES options for the cycle on level i
95640244768SBarry Smith -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
95740244768SBarry Smith 
95840244768SBarry Smith Notes:
95940244768SBarry Smith    The organization of the FAS solver is slightly different from the organization of PCMG
96040244768SBarry Smith    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
96140244768SBarry Smith    The cycle SNES instance may be used for monitoring convergence on a particular level.
96240244768SBarry Smith 
96340244768SBarry Smith Level: beginner
96440244768SBarry Smith 
96540244768SBarry Smith    References:
96640244768SBarry Smith . 1. -  Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
96740244768SBarry Smith    SIAM Review, 57(4), 2015
96840244768SBarry Smith 
96940244768SBarry Smith .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
97040244768SBarry Smith M*/
97140244768SBarry Smith 
97240244768SBarry Smith #undef __FUNCT__
97340244768SBarry Smith #define __FUNCT__ "SNESCreate_FAS"
97440244768SBarry Smith PETSC_EXTERN PetscErrorCode SNESCreate_FAS(SNES snes)
97540244768SBarry Smith {
97640244768SBarry Smith   SNES_FAS       *fas;
97740244768SBarry Smith   PetscErrorCode ierr;
97840244768SBarry Smith 
97940244768SBarry Smith   PetscFunctionBegin;
98040244768SBarry Smith   snes->ops->destroy        = SNESDestroy_FAS;
98140244768SBarry Smith   snes->ops->setup          = SNESSetUp_FAS;
98240244768SBarry Smith   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
98340244768SBarry Smith   snes->ops->view           = SNESView_FAS;
98440244768SBarry Smith   snes->ops->solve          = SNESSolve_FAS;
98540244768SBarry Smith   snes->ops->reset          = SNESReset_FAS;
98640244768SBarry Smith 
98740244768SBarry Smith   snes->usesksp = PETSC_FALSE;
98840244768SBarry Smith   snes->usespc  = PETSC_FALSE;
98940244768SBarry Smith 
99040244768SBarry Smith   if (!snes->tolerancesset) {
99140244768SBarry Smith     snes->max_funcs = 30000;
99240244768SBarry Smith     snes->max_its   = 10000;
99340244768SBarry Smith   }
99440244768SBarry Smith 
9954fc747eaSLawrence Mitchell   snes->alwayscomputesfinalresidual = PETSC_TRUE;
9964fc747eaSLawrence Mitchell 
99740244768SBarry Smith   ierr = PetscNewLog(snes,&fas);CHKERRQ(ierr);
99840244768SBarry Smith 
99940244768SBarry Smith   snes->data                  = (void*) fas;
100040244768SBarry Smith   fas->level                  = 0;
100140244768SBarry Smith   fas->levels                 = 1;
100240244768SBarry Smith   fas->n_cycles               = 1;
100340244768SBarry Smith   fas->max_up_it              = 1;
100440244768SBarry Smith   fas->max_down_it            = 1;
100540244768SBarry Smith   fas->smoothu                = NULL;
100640244768SBarry Smith   fas->smoothd                = NULL;
100740244768SBarry Smith   fas->next                   = NULL;
100840244768SBarry Smith   fas->previous               = NULL;
100940244768SBarry Smith   fas->fine                   = snes;
101040244768SBarry Smith   fas->interpolate            = NULL;
101140244768SBarry Smith   fas->restrct                = NULL;
101240244768SBarry Smith   fas->inject                 = NULL;
101340244768SBarry Smith   fas->usedmfornumberoflevels = PETSC_FALSE;
101440244768SBarry Smith   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
101540244768SBarry Smith   fas->full_downsweep         = PETSC_FALSE;
101640244768SBarry Smith 
101740244768SBarry Smith   fas->eventsmoothsetup    = 0;
101840244768SBarry Smith   fas->eventsmoothsolve    = 0;
101940244768SBarry Smith   fas->eventresidual       = 0;
102040244768SBarry Smith   fas->eventinterprestrict = 0;
102140244768SBarry Smith   PetscFunctionReturn(0);
102240244768SBarry Smith }
1023