xref: /petsc/src/snes/impls/fas/fas.c (revision b9c2fdf168d07ed48fb34d58b59f23d0c0351c9e)
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) {
173fde0ff24SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
17407144faaSPeter Brune   } else {
175ddebd997SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);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;
243162d76ddSPeter Brune   PetscBool      flg, upflg, downflg, monflg, galerkinflg;
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__
36439bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth"
36539bd7f45SPeter Brune /*
36639bd7f45SPeter Brune Defines the action of the downsmoother
36739bd7f45SPeter Brune  */
368*b9c2fdf1SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
369*b9c2fdf1SPeter 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);
376ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
377742fe5e2SPeter Brune   /* check convergence reason for the smoother */
378ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
379e70c42e5SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) {
380742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
381742fe5e2SPeter Brune     PetscFunctionReturn(0);
382742fe5e2SPeter Brune   }
383ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
3844b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
385*b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr);
38639bd7f45SPeter Brune   PetscFunctionReturn(0);
38739bd7f45SPeter Brune }
38839bd7f45SPeter Brune 
38939bd7f45SPeter Brune 
39039bd7f45SPeter Brune #undef __FUNCT__
39139bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth"
39239bd7f45SPeter Brune /*
39307144faaSPeter Brune Defines the action of the upsmoother
39439bd7f45SPeter Brune  */
395*b9c2fdf1SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) {
39639bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
39739bd7f45SPeter Brune   SNESConvergedReason reason;
398ab8d36c9SPeter Brune   Vec                 FPC;
399ab8d36c9SPeter Brune   SNES                smoothu;
40039bd7f45SPeter Brune   PetscFunctionBegin;
401ab8d36c9SPeter Brune 
402ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
403ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
40439bd7f45SPeter Brune   /* check convergence reason for the smoother */
405ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
40639bd7f45SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
40739bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
40839bd7f45SPeter Brune     PetscFunctionReturn(0);
40939bd7f45SPeter Brune   }
410ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
4114b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
412*b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr);
41339bd7f45SPeter Brune   PetscFunctionReturn(0);
41439bd7f45SPeter Brune }
41539bd7f45SPeter Brune 
41639bd7f45SPeter Brune #undef __FUNCT__
417938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
418938e4a01SJed Brown /*@
419938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
420938e4a01SJed Brown 
421938e4a01SJed Brown    Collective
422938e4a01SJed Brown 
423938e4a01SJed Brown    Input Arguments:
424938e4a01SJed Brown .  snes - SNESFAS
425938e4a01SJed Brown 
426938e4a01SJed Brown    Output Arguments:
427938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
428938e4a01SJed Brown 
429938e4a01SJed Brown    Level: developer
430938e4a01SJed Brown 
431938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
432938e4a01SJed Brown @*/
433938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
434938e4a01SJed Brown {
435938e4a01SJed Brown   PetscErrorCode ierr;
436938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
437938e4a01SJed Brown 
438938e4a01SJed Brown   PetscFunctionBegin;
439938e4a01SJed Brown   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
440938e4a01SJed Brown   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
441938e4a01SJed Brown   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
442938e4a01SJed Brown   PetscFunctionReturn(0);
443938e4a01SJed Brown }
444938e4a01SJed Brown 
445e9923e8dSJed Brown #undef __FUNCT__
446e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
447e9923e8dSJed Brown /*@
448e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
449e9923e8dSJed Brown 
450e9923e8dSJed Brown    Collective
451e9923e8dSJed Brown 
452e9923e8dSJed Brown    Input Arguments:
453e9923e8dSJed Brown +  fine - SNES from which to restrict
454e9923e8dSJed Brown -  Xfine - vector to restrict
455e9923e8dSJed Brown 
456e9923e8dSJed Brown    Output Arguments:
457e9923e8dSJed Brown .  Xcoarse - result of restriction
458e9923e8dSJed Brown 
459e9923e8dSJed Brown    Level: developer
460e9923e8dSJed Brown 
461e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
462e9923e8dSJed Brown @*/
463e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
464e9923e8dSJed Brown {
465e9923e8dSJed Brown   PetscErrorCode ierr;
466e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
467e9923e8dSJed Brown 
468e9923e8dSJed Brown   PetscFunctionBegin;
469e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
470e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
471e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
472e9923e8dSJed Brown   if (fas->inject) {
473e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
474e9923e8dSJed Brown   } else {
475e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
476e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
477e9923e8dSJed Brown   }
478e9923e8dSJed Brown   PetscFunctionReturn(0);
479e9923e8dSJed Brown }
480e9923e8dSJed Brown 
481e9923e8dSJed Brown #undef __FUNCT__
48239bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection"
48339bd7f45SPeter Brune /*
48439bd7f45SPeter Brune 
48539bd7f45SPeter Brune Performs the FAS coarse correction as:
48639bd7f45SPeter Brune 
48739bd7f45SPeter Brune fine problem: F(x) = 0
48839bd7f45SPeter Brune coarse problem: F^c(x) = b^c
48939bd7f45SPeter Brune 
49039bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
49139bd7f45SPeter Brune 
49239bd7f45SPeter Brune  */
49339a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
49439bd7f45SPeter Brune   PetscErrorCode      ierr;
49539bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
49639bd7f45SPeter Brune   SNESConvergedReason reason;
497ab8d36c9SPeter Brune   SNES                next;
498ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
49939bd7f45SPeter Brune   PetscFunctionBegin;
500ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
501ab8d36c9SPeter Brune   if (next) {
502ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
503ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
504ab8d36c9SPeter Brune 
505ab8d36c9SPeter Brune     X_c  = next->vec_sol;
506ab8d36c9SPeter Brune     Xo_c = next->work[0];
507ab8d36c9SPeter Brune     F_c  = next->vec_func;
508ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
509efe1f98aSPeter Brune 
510938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
511293a7e31SPeter Brune     /* restrict the defect */
512ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
513293a7e31SPeter Brune 
514*b9c2fdf1SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
515ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
516*b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
517c90fad12SPeter Brune 
518ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
519ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
520c90fad12SPeter Brune 
521c90fad12SPeter Brune     /* recurse to the next level */
522ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
523ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
524742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
525742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
526742fe5e2SPeter Brune       PetscFunctionReturn(0);
527742fe5e2SPeter Brune     }
528fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
529fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
530ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
531293a7e31SPeter Brune   }
53239bd7f45SPeter Brune   PetscFunctionReturn(0);
53339bd7f45SPeter Brune }
53439bd7f45SPeter Brune 
53539bd7f45SPeter Brune #undef __FUNCT__
53639bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive"
53739bd7f45SPeter Brune /*
53839bd7f45SPeter Brune 
53939bd7f45SPeter Brune The additive cycle looks like:
54039bd7f45SPeter Brune 
54107144faaSPeter Brune xhat = x
54207144faaSPeter Brune xhat = dS(x, b)
54307144faaSPeter Brune x = coarsecorrection(xhat, b_d)
54407144faaSPeter Brune x = x + nu*(xhat - x);
54539bd7f45SPeter Brune (optional) x = uS(x, b)
54639bd7f45SPeter Brune 
54739bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
54839bd7f45SPeter Brune 
54939bd7f45SPeter Brune  */
55039bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
55107144faaSPeter Brune   Vec                 F, B, Xhat;
55222c1e704SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
55339bd7f45SPeter Brune   PetscErrorCode      ierr;
55407144faaSPeter Brune   SNESConvergedReason reason;
55522c1e704SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
55622c1e704SPeter Brune   PetscBool           lssuccess;
557ab8d36c9SPeter Brune   SNES                next;
558ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
55939bd7f45SPeter Brune   PetscFunctionBegin;
560ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
56139bd7f45SPeter Brune   F = snes->vec_func;
56239bd7f45SPeter Brune   B = snes->vec_rhs;
563fde0ff24SPeter Brune   Xhat = snes->work[3];
56407144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
56507144faaSPeter Brune   /* recurse first */
566ab8d36c9SPeter Brune   if (next) {
567ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
568ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
56907144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
57039bd7f45SPeter Brune 
571ab8d36c9SPeter Brune     X_c  = next->vec_sol;
572ab8d36c9SPeter Brune     Xo_c = next->work[0];
573ab8d36c9SPeter Brune     F_c  = next->vec_func;
574ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
57539bd7f45SPeter Brune 
576938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
57707144faaSPeter Brune     /* restrict the defect */
578ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
57907144faaSPeter Brune 
58007144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
581ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
582*b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
58307144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
58407144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
58507144faaSPeter Brune 
58607144faaSPeter Brune     /* recurse */
587ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
58807144faaSPeter Brune 
58907144faaSPeter Brune     /* smooth on this level */
590*b9c2fdf1SPeter Brune     ierr = FASDownSmooth(snes, B, X, F, &fnorm);CHKERRQ(ierr);
59107144faaSPeter Brune 
592ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
59307144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
59407144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
59507144faaSPeter Brune       PetscFunctionReturn(0);
59607144faaSPeter Brune     }
59707144faaSPeter Brune 
59807144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
599c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
600ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
60107144faaSPeter Brune 
602ddebd997SPeter Brune     /* additive correction of the coarse direction*/
603f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
604f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
6059e764e56SPeter Brune     if (!lssuccess) {
6069e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6079e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6089e764e56SPeter Brune         PetscFunctionReturn(0);
6099e764e56SPeter Brune       }
6109e764e56SPeter Brune     }
611*b9c2fdf1SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
61207144faaSPeter Brune   } else {
613*b9c2fdf1SPeter Brune     ierr = FASDownSmooth(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
61407144faaSPeter Brune   }
61539bd7f45SPeter Brune   PetscFunctionReturn(0);
61639bd7f45SPeter Brune }
61739bd7f45SPeter Brune 
61839bd7f45SPeter Brune #undef __FUNCT__
61939bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative"
62039bd7f45SPeter Brune /*
62139bd7f45SPeter Brune 
62239bd7f45SPeter Brune Defines the FAS cycle as:
62339bd7f45SPeter Brune 
62439bd7f45SPeter Brune fine problem: F(x) = 0
62539bd7f45SPeter Brune coarse problem: F^c(x) = b^c
62639bd7f45SPeter Brune 
62739bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
62839bd7f45SPeter Brune 
62939bd7f45SPeter Brune correction:
63039bd7f45SPeter Brune 
63139bd7f45SPeter Brune x = x + I(x^c - Rx)
63239bd7f45SPeter Brune 
63339bd7f45SPeter Brune  */
63439bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
63539bd7f45SPeter Brune 
63639bd7f45SPeter Brune   PetscErrorCode      ierr;
63739bd7f45SPeter Brune   Vec                 F,B;
63839bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
63939bd7f45SPeter Brune 
64039bd7f45SPeter Brune   PetscFunctionBegin;
64139bd7f45SPeter Brune   F = snes->vec_func;
64239bd7f45SPeter Brune   B = snes->vec_rhs;
64339bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
644*b9c2fdf1SPeter Brune   ierr = FASDownSmooth(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
64539bd7f45SPeter Brune 
646c90fad12SPeter Brune   if (fas->level != 0) {
647e70c42e5SPeter Brune     ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
648*b9c2fdf1SPeter Brune     ierr = FASUpSmooth(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
649fe6f9142SPeter Brune   }
650fa9694d7SPeter Brune   PetscFunctionReturn(0);
651421d9b32SPeter Brune }
652421d9b32SPeter Brune 
653421d9b32SPeter Brune #undef __FUNCT__
654421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
655421d9b32SPeter Brune 
656421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
657421d9b32SPeter Brune {
658fa9694d7SPeter Brune   PetscErrorCode ierr;
659fe6f9142SPeter Brune   PetscInt       i, maxits;
660ddb5aff1SPeter Brune   Vec            X, F;
661fe6f9142SPeter Brune   PetscReal      fnorm;
662b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
663b17ce1afSJed Brown   DM             dm;
664e70c42e5SPeter Brune   PetscBool      isFine;
665b17ce1afSJed Brown 
666421d9b32SPeter Brune   PetscFunctionBegin;
667fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
668fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
669fa9694d7SPeter Brune   X = snes->vec_sol;
670f5a6d4f9SBarry Smith   F = snes->vec_func;
671293a7e31SPeter Brune 
67218a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
673293a7e31SPeter Brune   /*norm setup */
674fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
675fe6f9142SPeter Brune   snes->iter = 0;
676fe6f9142SPeter Brune   snes->norm = 0.;
677fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
67818a66777SPeter Brune   if (isFine || fas->monitor) {
679fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
680fe6f9142SPeter Brune     if (snes->domainerror) {
681fe6f9142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
682fe6f9142SPeter Brune       PetscFunctionReturn(0);
683fe6f9142SPeter Brune     }
684fe6f9142SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
685fe6f9142SPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
686fe6f9142SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
687fe6f9142SPeter Brune     snes->norm = fnorm;
688fe6f9142SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
689fe6f9142SPeter Brune     SNESLogConvHistory(snes,fnorm,0);
690fe6f9142SPeter Brune     ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
691fe6f9142SPeter Brune 
692fe6f9142SPeter Brune     /* set parameter for default relative tolerance convergence test */
693fe6f9142SPeter Brune     snes->ttol = fnorm*snes->rtol;
694fe6f9142SPeter Brune     /* test convergence */
695fe6f9142SPeter Brune     ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
696fe6f9142SPeter Brune     if (snes->reason) PetscFunctionReturn(0);
69766585501SPeter Brune   }
698b17ce1afSJed Brown 
699*b9c2fdf1SPeter Brune   if (isFine) {
700*b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
701b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
702b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
703b17ce1afSJed Brown       DM dmcoarse;
704b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
705b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
706b17ce1afSJed Brown       dm = dmcoarse;
707b17ce1afSJed Brown     }
708*b9c2fdf1SPeter Brune   }
709b17ce1afSJed Brown 
710fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
711fe6f9142SPeter Brune     /* Call general purpose update function */
712646217ecSPeter Brune 
713fe6f9142SPeter Brune     if (snes->ops->update) {
714fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
715fe6f9142SPeter Brune     }
71607144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
71739bd7f45SPeter Brune       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
71807144faaSPeter Brune     } else {
71907144faaSPeter Brune       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
72007144faaSPeter Brune     }
721742fe5e2SPeter Brune 
722742fe5e2SPeter Brune     /* check for FAS cycle divergence */
723742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
724742fe5e2SPeter Brune       PetscFunctionReturn(0);
725742fe5e2SPeter Brune     }
726*b9c2fdf1SPeter Brune 
727c90fad12SPeter Brune     /* Monitor convergence */
728c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
729c90fad12SPeter Brune     snes->iter = i+1;
730c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
731c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
732c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
733c90fad12SPeter Brune     /* Test for convergence */
73466585501SPeter Brune     if (isFine) {
735*b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
736c90fad12SPeter Brune       if (snes->reason) break;
737fe6f9142SPeter Brune     }
73866585501SPeter Brune   }
739fe6f9142SPeter Brune   if (i == maxits) {
740fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
741fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
742fe6f9142SPeter Brune   }
743421d9b32SPeter Brune   PetscFunctionReturn(0);
744421d9b32SPeter Brune }
745