xref: /petsc/src/snes/impls/fas/fas.c (revision e4ed7901070ce1ca36eb7d68dd07557542e66e31)
1421d9b32SPeter Brune /* Defines the basic SNES object */
26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3421d9b32SPeter Brune 
407144faaSPeter Brune const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0};
507144faaSPeter Brune 
6421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
7421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
8421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
9421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
10421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
11421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
126273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *);
13ab8d36c9SPeter Brune extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *);
14421d9b32SPeter Brune 
151fbfccc6SJed Brown /*MC
161fbfccc6SJed Brown 
171fbfccc6SJed Brown SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
181fbfccc6SJed Brown 
19d3bc2379SPeter Brune    The nonlinear problem is solved by correction using coarse versions
20d3bc2379SPeter Brune    of the nonlinear problem.  This problem is perturbed so that a projected
21d3bc2379SPeter Brune    solution of the fine problem elicits no correction from the coarse problem.
22d3bc2379SPeter Brune 
23d3bc2379SPeter Brune Options Database:
24d3bc2379SPeter Brune +   -snes_fas_levels -  The number of levels
25d3bc2379SPeter Brune .   -snes_fas_cycles<1> -  The number of cycles -- 1 for V, 2 for W
26d3bc2379SPeter Brune .   -snes_fas_type<additive, multiplicative>  -  Additive or multiplicative cycle
27d3bc2379SPeter Brune .   -snes_fas_galerkin<PETSC_FALSE> -  Form coarse problems by projection back upon the fine problem
28d3bc2379SPeter Brune .   -snes_fas_smoothup<1> -  The number of iterations of the post-smoother
29d3bc2379SPeter Brune .   -snes_fas_smoothdown<1> -  The number of iterations of the pre-smoother
30d3bc2379SPeter Brune .   -snes_fas_monitor -  Monitor progress of all of the levels
31d3bc2379SPeter Brune .   -fas_levels_snes_ -  SNES options for all smoothers
32d3bc2379SPeter Brune .   -fas_levels_cycle_snes -  SNES options for all cycles
33d3bc2379SPeter Brune .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
34d3bc2379SPeter Brune .   -fas_levels_i_cycle_snes - SNES options for the cycle on level i
35d3bc2379SPeter Brune -   -fas_coarse_snes_ -  SNES options for the coarsest smoother
36d3bc2379SPeter Brune 
37d3bc2379SPeter Brune Notes:
38d3bc2379SPeter Brune    The organization of the FAS solver is slightly different from the organization of PCMG
39d3bc2379SPeter Brune    As each level has smoother SNES instances(down and potentially up) and a cycle SNES instance.
40d3bc2379SPeter Brune    The cycle SNES instance may be used for monitoring convergence on a particular level.
411fbfccc6SJed Brown 
421fbfccc6SJed Brown Level: advanced
431fbfccc6SJed Brown 
44d3bc2379SPeter Brune .seealso: PCMG, SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
451fbfccc6SJed Brown M*/
46421d9b32SPeter Brune 
47421d9b32SPeter Brune #undef __FUNCT__
48421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
491fbfccc6SJed Brown PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes)
50421d9b32SPeter Brune {
51421d9b32SPeter Brune   SNES_FAS * fas;
52421d9b32SPeter Brune   PetscErrorCode ierr;
53421d9b32SPeter Brune 
54421d9b32SPeter Brune   PetscFunctionBegin;
55421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
56421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
57421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
58421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
59421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
60421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
61421d9b32SPeter Brune 
62ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
63ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
64ed020824SBarry Smith 
6588976e71SPeter Brune   if (!snes->tolerancesset) {
660e444f03SPeter Brune     snes->max_funcs = 30000;
670e444f03SPeter Brune     snes->max_its   = 10000;
6888976e71SPeter Brune   }
690e444f03SPeter Brune 
70421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
71421d9b32SPeter Brune   snes->data                  = (void*) fas;
72421d9b32SPeter Brune   fas->level                  = 0;
73293a7e31SPeter Brune   fas->levels                 = 1;
74ee78dd50SPeter Brune   fas->n_cycles               = 1;
75ee78dd50SPeter Brune   fas->max_up_it              = 1;
76ee78dd50SPeter Brune   fas->max_down_it            = 1;
77ab8d36c9SPeter Brune   fas->smoothu                = PETSC_NULL;
78ab8d36c9SPeter Brune   fas->smoothd                = PETSC_NULL;
79421d9b32SPeter Brune   fas->next                   = PETSC_NULL;
806273346dSPeter Brune   fas->previous               = PETSC_NULL;
81ab8d36c9SPeter Brune   fas->fine                   = snes;
82421d9b32SPeter Brune   fas->interpolate            = PETSC_NULL;
83421d9b32SPeter Brune   fas->restrct                = PETSC_NULL;
84efe1f98aSPeter Brune   fas->inject                 = PETSC_NULL;
85cc05f883SPeter Brune   fas->monitor                = PETSC_NULL;
86cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
87ddebd997SPeter Brune   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
88efe1f98aSPeter Brune   PetscFunctionReturn(0);
89efe1f98aSPeter Brune }
90efe1f98aSPeter Brune 
91421d9b32SPeter Brune #undef __FUNCT__
92421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
93421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
94421d9b32SPeter Brune {
9577df8cc4SPeter Brune   PetscErrorCode ierr = 0;
96421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
97421d9b32SPeter Brune 
98421d9b32SPeter Brune   PetscFunctionBegin;
99ab8d36c9SPeter Brune   if (fas->smoothu != fas->smoothd) {
100ab8d36c9SPeter Brune     ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
101ab8d36c9SPeter Brune     ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
102ab8d36c9SPeter Brune   } else {
103ab8d36c9SPeter Brune     ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
104ab8d36c9SPeter Brune     fas->smoothu = PETSC_NULL;
105ab8d36c9SPeter Brune   }
1063dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
107bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
108bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
109bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
110742fe5e2SPeter Brune   if (fas->next)      ierr = SNESReset(fas->next);CHKERRQ(ierr);
11122c1e704SPeter Brune 
112421d9b32SPeter Brune   PetscFunctionReturn(0);
113421d9b32SPeter Brune }
114421d9b32SPeter Brune 
115421d9b32SPeter Brune #undef __FUNCT__
116421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
117421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
118421d9b32SPeter Brune {
119421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
120742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
121421d9b32SPeter Brune 
122421d9b32SPeter Brune   PetscFunctionBegin;
123421d9b32SPeter Brune   /* recursively resets and then destroys */
12479d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
125421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
126421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
127ee1fd11aSPeter Brune 
128421d9b32SPeter Brune   PetscFunctionReturn(0);
129421d9b32SPeter Brune }
130421d9b32SPeter Brune 
131421d9b32SPeter Brune #undef __FUNCT__
132421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
133421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
134421d9b32SPeter Brune {
13548bfdf8aSPeter Brune   SNES_FAS                *fas = (SNES_FAS *) snes->data;
136421d9b32SPeter Brune   PetscErrorCode          ierr;
137efe1f98aSPeter Brune   VecScatter              injscatter;
138d1adcc6fSPeter Brune   PetscInt                dm_levels;
1393dccd265SPeter Brune   Vec                     vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
140ab8d36c9SPeter Brune   SNES                    next;
141ab8d36c9SPeter Brune   PetscBool               isFine;
142421d9b32SPeter Brune   PetscFunctionBegin;
143eff52c0eSPeter Brune 
144ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
145ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
146d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
147d1adcc6fSPeter Brune     dm_levels++;
148cc05f883SPeter Brune     if (dm_levels > fas->levels) {
1492e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
1503dccd265SPeter Brune       vec_sol = snes->vec_sol;
1513dccd265SPeter Brune       vec_func = snes->vec_func;
1523dccd265SPeter Brune       vec_sol_update = snes->vec_sol_update;
1533dccd265SPeter Brune       vec_rhs = snes->vec_rhs;
1543dccd265SPeter Brune       snes->vec_sol = PETSC_NULL;
1553dccd265SPeter Brune       snes->vec_func = PETSC_NULL;
1563dccd265SPeter Brune       snes->vec_sol_update = PETSC_NULL;
1573dccd265SPeter Brune       snes->vec_rhs = PETSC_NULL;
1583dccd265SPeter Brune 
1593dccd265SPeter Brune       /* reset the number of levels */
160d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
161cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
1623dccd265SPeter Brune 
1633dccd265SPeter Brune       snes->vec_sol = vec_sol;
1643dccd265SPeter Brune       snes->vec_func = vec_func;
1653dccd265SPeter Brune       snes->vec_rhs = vec_rhs;
1663dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
167d1adcc6fSPeter Brune     }
168d1adcc6fSPeter Brune   }
169ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
170ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
1713dccd265SPeter Brune 
17207144faaSPeter Brune   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
173*e4ed7901SPeter Brune     ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
17407144faaSPeter Brune   } else {
175e7f468e7SPeter Brune     ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
17607144faaSPeter Brune   }
177cc05f883SPeter Brune 
178ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
179ab8d36c9SPeter Brune   if (!fas->smoothd) {
180ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
181ab8d36c9SPeter Brune   }
182ab8d36c9SPeter Brune 
18379d9a41aSPeter Brune   if (snes->dm) {
184ab8d36c9SPeter Brune     /* set the smoother DMs properly */
185ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
186ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
18779d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
188ab8d36c9SPeter Brune     if (next) {
18979d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
190ab8d36c9SPeter Brune       if (!next->dm) {
191ab8d36c9SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr);
192ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
19379d9a41aSPeter Brune       }
19479d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
19579d9a41aSPeter Brune       if (!fas->interpolate) {
196ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
197bccf9bb3SJed Brown         if (!fas->restrct) {
198bccf9bb3SJed Brown           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
19979d9a41aSPeter Brune           fas->restrct = fas->interpolate;
20079d9a41aSPeter Brune         }
201bccf9bb3SJed Brown       }
20279d9a41aSPeter Brune       /* set the injection from the DM */
20379d9a41aSPeter Brune       if (!fas->inject) {
204ab8d36c9SPeter Brune         ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
20579d9a41aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
20679d9a41aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
20779d9a41aSPeter Brune       }
20879d9a41aSPeter Brune     }
20979d9a41aSPeter Brune   }
21079d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
21179d9a41aSPeter Brune 
21279d9a41aSPeter Brune   if (fas->galerkin) {
213ab8d36c9SPeter Brune     if (next)
214ab8d36c9SPeter Brune       ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
215ab8d36c9SPeter Brune     if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
216ab8d36c9SPeter Brune     if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
21779d9a41aSPeter Brune   }
21879d9a41aSPeter Brune 
219d06165b7SPeter Brune   /* set the smoothers up here so that precedence is taken for instance-specific options over the whole-solver options */
220d06165b7SPeter Brune   if(fas->smoothu) ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
221d06165b7SPeter Brune   if(fas->smoothd) ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
222d06165b7SPeter Brune 
223ab8d36c9SPeter Brune   if (next) {
22479d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
225ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
226ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
227ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
22879d9a41aSPeter Brune   }
2296273346dSPeter Brune   /* setup FAS work vectors */
2306273346dSPeter Brune   if (fas->galerkin) {
2316273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
2326273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
2336273346dSPeter Brune   }
234421d9b32SPeter Brune   PetscFunctionReturn(0);
235421d9b32SPeter Brune }
236421d9b32SPeter Brune 
237421d9b32SPeter Brune #undef __FUNCT__
238421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
239421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
240421d9b32SPeter Brune {
241ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
242ee78dd50SPeter Brune   PetscInt       levels = 1;
2434d26bfa5SPeter Brune   PetscBool      flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE;
244421d9b32SPeter Brune   PetscErrorCode ierr;
245ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
24607144faaSPeter Brune   SNESFASType    fastype;
247fde0ff24SPeter Brune   const char     *optionsprefix;
248f1c6b773SPeter Brune   SNESLineSearch linesearch;
24966585501SPeter Brune   PetscInt       m, n_up, n_down;
250ab8d36c9SPeter Brune   SNES           next;
251ab8d36c9SPeter Brune   PetscBool      isFine;
252421d9b32SPeter Brune 
253421d9b32SPeter Brune   PetscFunctionBegin;
254ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
255c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
256ee78dd50SPeter Brune 
257ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
258ab8d36c9SPeter Brune   if (isFine) {
259ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
260c732cbdbSBarry Smith     if (!flg && snes->dm) {
261c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
262c732cbdbSBarry Smith       levels++;
263d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
264c732cbdbSBarry Smith     }
265ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
26607144faaSPeter Brune     fastype = fas->fastype;
26707144faaSPeter Brune     ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
26807144faaSPeter Brune     if (flg) {
26907144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
27007144faaSPeter Brune     }
271ee78dd50SPeter Brune 
272fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
273ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
274ab8d36c9SPeter Brune     if (flg) {
275ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
276fde0ff24SPeter Brune     }
277fde0ff24SPeter Brune 
278ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
279ab8d36c9SPeter Brune     if (flg) {
280ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
281ab8d36c9SPeter Brune     }
282ee78dd50SPeter Brune 
28366585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
284162d76ddSPeter Brune 
28566585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
286162d76ddSPeter Brune 
287c8c899caSPeter Brune     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
288c8c899caSPeter Brune     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
289ab8d36c9SPeter Brune   }
290ee78dd50SPeter Brune 
291421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
2928cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
293162d76ddSPeter Brune   if (upflg) {
29466585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
295162d76ddSPeter Brune   }
296162d76ddSPeter Brune   if (downflg) {
29766585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
298162d76ddSPeter Brune   }
299eff52c0eSPeter Brune 
3009e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
3019e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
3029e764e56SPeter Brune     if (!snes->linesearch) {
303f1c6b773SPeter Brune       ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
304f1c6b773SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr);
3059e764e56SPeter Brune     }
3069e764e56SPeter Brune   }
3079e764e56SPeter Brune 
308ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
309ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
310ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
311421d9b32SPeter Brune   PetscFunctionReturn(0);
312421d9b32SPeter Brune }
313421d9b32SPeter Brune 
314421d9b32SPeter Brune #undef __FUNCT__
315421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
316421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
317421d9b32SPeter Brune {
318421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
319ab8d36c9SPeter Brune   PetscBool      isFine, iascii;
320ab8d36c9SPeter Brune   PetscInt       i;
321421d9b32SPeter Brune   PetscErrorCode ierr;
322ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
323421d9b32SPeter Brune 
324421d9b32SPeter Brune   PetscFunctionBegin;
325ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
326ab8d36c9SPeter Brune   if (isFine) {
327421d9b32SPeter Brune     ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
328421d9b32SPeter Brune     if (iascii) {
329ab8d36c9SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
330ab8d36c9SPeter Brune       if (fas->galerkin) {
331ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
332421d9b32SPeter Brune       } else {
333ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
334421d9b32SPeter Brune       }
335ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
336ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
337ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
338ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
339ab8d36c9SPeter Brune         if (!i) {
340ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
341421d9b32SPeter Brune         } else {
342ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
343421d9b32SPeter Brune         }
344ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
345ab8d36c9SPeter Brune         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
346ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
347ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
348ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
349ab8d36c9SPeter Brune         } else if (i) {
350ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
351ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
352ab8d36c9SPeter Brune           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
353ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
354ab8d36c9SPeter Brune         }
355ab8d36c9SPeter Brune       }
356421d9b32SPeter Brune     } else {
357421d9b32SPeter Brune       SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
358421d9b32SPeter Brune     }
359ab8d36c9SPeter Brune   }
360421d9b32SPeter Brune   PetscFunctionReturn(0);
361421d9b32SPeter Brune }
362421d9b32SPeter Brune 
363421d9b32SPeter Brune #undef __FUNCT__
36491f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private"
36539bd7f45SPeter Brune /*
36639bd7f45SPeter Brune Defines the action of the downsmoother
36739bd7f45SPeter Brune  */
36891f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
369b9c2fdf1SPeter Brune {
37039bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
371742fe5e2SPeter Brune   SNESConvergedReason reason;
372ab8d36c9SPeter Brune   Vec                 FPC;
373ab8d36c9SPeter Brune   SNES                smoothd;
374421d9b32SPeter Brune   PetscFunctionBegin;
375ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
376*e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
377*e4ed7901SPeter Brune   ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr);
378ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
379742fe5e2SPeter Brune   /* check convergence reason for the smoother */
380ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
381e70c42e5SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) {
382742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
383742fe5e2SPeter Brune     PetscFunctionReturn(0);
384742fe5e2SPeter Brune   }
385ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
3864b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
387b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr);
38839bd7f45SPeter Brune   PetscFunctionReturn(0);
38939bd7f45SPeter Brune }
39039bd7f45SPeter Brune 
39139bd7f45SPeter Brune 
39239bd7f45SPeter Brune #undef __FUNCT__
39391f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private"
39439bd7f45SPeter Brune /*
39507144faaSPeter Brune Defines the action of the upsmoother
39639bd7f45SPeter Brune  */
39791f99d7cSPeter Brune PetscErrorCode SNESFASUpSmooth_Private (SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) {
39839bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
39939bd7f45SPeter Brune   SNESConvergedReason reason;
400ab8d36c9SPeter Brune   Vec                 FPC;
401ab8d36c9SPeter Brune   SNES                smoothu;
40239bd7f45SPeter Brune   PetscFunctionBegin;
403ab8d36c9SPeter Brune 
404ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
405ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
40639bd7f45SPeter Brune   /* check convergence reason for the smoother */
407ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
40839bd7f45SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
40939bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
41039bd7f45SPeter Brune     PetscFunctionReturn(0);
41139bd7f45SPeter Brune   }
412ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
4134b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
414b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr);
41539bd7f45SPeter Brune   PetscFunctionReturn(0);
41639bd7f45SPeter Brune }
41739bd7f45SPeter Brune 
41839bd7f45SPeter Brune #undef __FUNCT__
419938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
420938e4a01SJed Brown /*@
421938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
422938e4a01SJed Brown 
423938e4a01SJed Brown    Collective
424938e4a01SJed Brown 
425938e4a01SJed Brown    Input Arguments:
426938e4a01SJed Brown .  snes - SNESFAS
427938e4a01SJed Brown 
428938e4a01SJed Brown    Output Arguments:
429938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
430938e4a01SJed Brown 
431938e4a01SJed Brown    Level: developer
432938e4a01SJed Brown 
433938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
434938e4a01SJed Brown @*/
435938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
436938e4a01SJed Brown {
437938e4a01SJed Brown   PetscErrorCode ierr;
438938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
439938e4a01SJed Brown 
440938e4a01SJed Brown   PetscFunctionBegin;
441938e4a01SJed Brown   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
442938e4a01SJed Brown   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
443938e4a01SJed Brown   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
444938e4a01SJed Brown   PetscFunctionReturn(0);
445938e4a01SJed Brown }
446938e4a01SJed Brown 
447e9923e8dSJed Brown #undef __FUNCT__
448e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
449e9923e8dSJed Brown /*@
450e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
451e9923e8dSJed Brown 
452e9923e8dSJed Brown    Collective
453e9923e8dSJed Brown 
454e9923e8dSJed Brown    Input Arguments:
455e9923e8dSJed Brown +  fine - SNES from which to restrict
456e9923e8dSJed Brown -  Xfine - vector to restrict
457e9923e8dSJed Brown 
458e9923e8dSJed Brown    Output Arguments:
459e9923e8dSJed Brown .  Xcoarse - result of restriction
460e9923e8dSJed Brown 
461e9923e8dSJed Brown    Level: developer
462e9923e8dSJed Brown 
463e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
464e9923e8dSJed Brown @*/
465e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
466e9923e8dSJed Brown {
467e9923e8dSJed Brown   PetscErrorCode ierr;
468e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
469e9923e8dSJed Brown 
470e9923e8dSJed Brown   PetscFunctionBegin;
471e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
472e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
473e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
474e9923e8dSJed Brown   if (fas->inject) {
475e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
476e9923e8dSJed Brown   } else {
477e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
478e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
479e9923e8dSJed Brown   }
480e9923e8dSJed Brown   PetscFunctionReturn(0);
481e9923e8dSJed Brown }
482e9923e8dSJed Brown 
483e9923e8dSJed Brown #undef __FUNCT__
4848c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection"
48539bd7f45SPeter Brune /*
48639bd7f45SPeter Brune 
48739bd7f45SPeter Brune Performs the FAS coarse correction as:
48839bd7f45SPeter Brune 
48939bd7f45SPeter Brune fine problem: F(x) = 0
49039bd7f45SPeter Brune coarse problem: F^c(x) = b^c
49139bd7f45SPeter Brune 
49239bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
49339bd7f45SPeter Brune 
49439bd7f45SPeter Brune  */
4958c40d5fbSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
49639bd7f45SPeter Brune   PetscErrorCode      ierr;
49739bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
49839bd7f45SPeter Brune   SNESConvergedReason reason;
499ab8d36c9SPeter Brune   SNES                next;
500ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
50139bd7f45SPeter Brune   PetscFunctionBegin;
502ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
503ab8d36c9SPeter Brune   if (next) {
504ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
505ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
506ab8d36c9SPeter Brune 
507ab8d36c9SPeter Brune     X_c  = next->vec_sol;
508ab8d36c9SPeter Brune     Xo_c = next->work[0];
509ab8d36c9SPeter Brune     F_c  = next->vec_func;
510ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
511efe1f98aSPeter Brune 
512938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
513293a7e31SPeter Brune     /* restrict the defect */
514ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
515b9c2fdf1SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
516ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
517*e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
518b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
519*e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
520ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
521ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
522c90fad12SPeter Brune 
523c90fad12SPeter Brune     /* recurse to the next level */
524*e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
525ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
526ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
527742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
528742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
529742fe5e2SPeter Brune       PetscFunctionReturn(0);
530742fe5e2SPeter Brune     }
531fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
532fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
533ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
534293a7e31SPeter Brune   }
53539bd7f45SPeter Brune   PetscFunctionReturn(0);
53639bd7f45SPeter Brune }
53739bd7f45SPeter Brune 
53839bd7f45SPeter Brune #undef __FUNCT__
5392cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive"
54039bd7f45SPeter Brune /*
54139bd7f45SPeter Brune 
54239bd7f45SPeter Brune The additive cycle looks like:
54339bd7f45SPeter Brune 
54407144faaSPeter Brune xhat = x
54507144faaSPeter Brune xhat = dS(x, b)
54607144faaSPeter Brune x = coarsecorrection(xhat, b_d)
54707144faaSPeter Brune x = x + nu*(xhat - x);
54839bd7f45SPeter Brune (optional) x = uS(x, b)
54939bd7f45SPeter Brune 
55039bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
55139bd7f45SPeter Brune 
55239bd7f45SPeter Brune  */
55391f99d7cSPeter Brune PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) {
55407144faaSPeter Brune   Vec                 F, B, Xhat;
55522c1e704SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
55639bd7f45SPeter Brune   PetscErrorCode      ierr;
55707144faaSPeter Brune   SNESConvergedReason reason;
55822c1e704SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
55922c1e704SPeter Brune   PetscBool           lssuccess;
560ab8d36c9SPeter Brune   SNES                next;
561ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
56239bd7f45SPeter Brune   PetscFunctionBegin;
563ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
56439bd7f45SPeter Brune   F = snes->vec_func;
56539bd7f45SPeter Brune   B = snes->vec_rhs;
566e7f468e7SPeter Brune   Xhat = snes->work[1];
56707144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
56807144faaSPeter Brune   /* recurse first */
569ab8d36c9SPeter Brune   if (next) {
570ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
571ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
57207144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
573*e4ed7901SPeter Brune     ierr = VecNorm(F, &fnorm);CHKERRQ(ierr);
574ab8d36c9SPeter Brune     X_c  = next->vec_sol;
575ab8d36c9SPeter Brune     Xo_c = next->work[0];
576ab8d36c9SPeter Brune     F_c  = next->vec_func;
577ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
57839bd7f45SPeter Brune 
579938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
58007144faaSPeter Brune     /* restrict the defect */
581ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
58207144faaSPeter Brune 
58307144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
584ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
585*e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
586b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
587*e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
58807144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
58907144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
59007144faaSPeter Brune 
59107144faaSPeter Brune     /* recurse */
592*e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
593ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
59407144faaSPeter Brune 
59507144faaSPeter Brune     /* smooth on this level */
59691f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
59707144faaSPeter Brune 
598ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
59907144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
60007144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
60107144faaSPeter Brune       PetscFunctionReturn(0);
60207144faaSPeter Brune     }
60307144faaSPeter Brune 
60407144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
605c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
606ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
60707144faaSPeter Brune 
608ddebd997SPeter Brune     /* additive correction of the coarse direction*/
609f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
610f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
6119e764e56SPeter Brune     if (!lssuccess) {
6129e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6139e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6149e764e56SPeter Brune         PetscFunctionReturn(0);
6159e764e56SPeter Brune       }
6169e764e56SPeter Brune     }
617b9c2fdf1SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
61807144faaSPeter Brune   } else {
61991f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
62007144faaSPeter Brune   }
62139bd7f45SPeter Brune   PetscFunctionReturn(0);
62239bd7f45SPeter Brune }
62339bd7f45SPeter Brune 
62439bd7f45SPeter Brune #undef __FUNCT__
6252cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative"
62639bd7f45SPeter Brune /*
62739bd7f45SPeter Brune 
62839bd7f45SPeter Brune Defines the FAS cycle as:
62939bd7f45SPeter Brune 
63039bd7f45SPeter Brune fine problem: F(x) = 0
63139bd7f45SPeter Brune coarse problem: F^c(x) = b^c
63239bd7f45SPeter Brune 
63339bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
63439bd7f45SPeter Brune 
63539bd7f45SPeter Brune correction:
63639bd7f45SPeter Brune 
63739bd7f45SPeter Brune x = x + I(x^c - Rx)
63839bd7f45SPeter Brune 
63939bd7f45SPeter Brune  */
64091f99d7cSPeter Brune PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) {
64139bd7f45SPeter Brune 
64239bd7f45SPeter Brune   PetscErrorCode      ierr;
64339bd7f45SPeter Brune   Vec                 F,B;
64439bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
64539bd7f45SPeter Brune 
64639bd7f45SPeter Brune   PetscFunctionBegin;
64739bd7f45SPeter Brune   F = snes->vec_func;
64839bd7f45SPeter Brune   B = snes->vec_rhs;
64939bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
65091f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
651c90fad12SPeter Brune   if (fas->level != 0) {
6528c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
65391f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
654fe6f9142SPeter Brune   }
655fa9694d7SPeter Brune   PetscFunctionReturn(0);
656421d9b32SPeter Brune }
657421d9b32SPeter Brune 
658421d9b32SPeter Brune #undef __FUNCT__
659421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
660421d9b32SPeter Brune 
661421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
662421d9b32SPeter Brune {
663fa9694d7SPeter Brune   PetscErrorCode ierr;
664fe6f9142SPeter Brune   PetscInt       i, maxits;
665ddb5aff1SPeter Brune   Vec            X, F;
666fe6f9142SPeter Brune   PetscReal      fnorm;
667b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
668b17ce1afSJed Brown   DM             dm;
669e70c42e5SPeter Brune   PetscBool      isFine;
670b17ce1afSJed Brown 
671421d9b32SPeter Brune   PetscFunctionBegin;
672fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
673fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
674fa9694d7SPeter Brune   X = snes->vec_sol;
675f5a6d4f9SBarry Smith   F = snes->vec_func;
676293a7e31SPeter Brune 
67718a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
678293a7e31SPeter Brune   /*norm setup */
679fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
680fe6f9142SPeter Brune   snes->iter = 0;
681fe6f9142SPeter Brune   snes->norm = 0.;
682fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
683*e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
684fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
685fe6f9142SPeter Brune     if (snes->domainerror) {
686fe6f9142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
687fe6f9142SPeter Brune       PetscFunctionReturn(0);
688fe6f9142SPeter Brune     }
689*e4ed7901SPeter Brune   } else {
690*e4ed7901SPeter Brune     snes->vec_func_init_set = PETSC_FALSE;
691*e4ed7901SPeter Brune   }
692*e4ed7901SPeter Brune 
693*e4ed7901SPeter Brune   if (!snes->norm_init_set) {
694fe6f9142SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
695fe6f9142SPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
696fe6f9142SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
697*e4ed7901SPeter Brune   } else {
698*e4ed7901SPeter Brune     fnorm = snes->norm_init;
699*e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
700*e4ed7901SPeter Brune   }
701*e4ed7901SPeter Brune 
702fe6f9142SPeter Brune   snes->norm = fnorm;
703fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
704fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
705fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
706fe6f9142SPeter Brune 
707fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
708fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
709fe6f9142SPeter Brune   /* test convergence */
710fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
711fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
712*e4ed7901SPeter Brune 
713b17ce1afSJed Brown 
714b9c2fdf1SPeter Brune   if (isFine) {
715b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
716b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
717b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
718b17ce1afSJed Brown       DM dmcoarse;
719b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
720b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
721b17ce1afSJed Brown       dm = dmcoarse;
722b17ce1afSJed Brown     }
723b9c2fdf1SPeter Brune   }
724b17ce1afSJed Brown 
725fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
726fe6f9142SPeter Brune     /* Call general purpose update function */
727646217ecSPeter Brune 
728fe6f9142SPeter Brune     if (snes->ops->update) {
729fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
730fe6f9142SPeter Brune     }
73107144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
73291f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
73307144faaSPeter Brune     } else {
73491f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
73507144faaSPeter Brune     }
736742fe5e2SPeter Brune 
737742fe5e2SPeter Brune     /* check for FAS cycle divergence */
738742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
739742fe5e2SPeter Brune       PetscFunctionReturn(0);
740742fe5e2SPeter Brune     }
741b9c2fdf1SPeter Brune 
742c90fad12SPeter Brune     /* Monitor convergence */
743c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
744c90fad12SPeter Brune     snes->iter = i+1;
745c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
746c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
747c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
748c90fad12SPeter Brune     /* Test for convergence */
74966585501SPeter Brune     if (isFine) {
750b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
751c90fad12SPeter Brune       if (snes->reason) break;
752fe6f9142SPeter Brune     }
75366585501SPeter Brune   }
754fe6f9142SPeter Brune   if (i == maxits) {
755fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
756fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
757fe6f9142SPeter Brune   }
758421d9b32SPeter Brune   PetscFunctionReturn(0);
759421d9b32SPeter Brune }
760