xref: /petsc/src/snes/impls/fas/fas.c (revision f89ba88e4820770605f92f5e73984a1592b2e9a1)
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
327d84e935SPeter Brune .   -fas_levels_cycle_snes_ -  SNES options for all cycles
33d3bc2379SPeter Brune .   -fas_levels_i_snes_ -  SNES options for the smoothers on level i
347d84e935SPeter 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 
427d84e935SPeter Brune Level: beginner
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   ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
100ab8d36c9SPeter Brune   ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
1013dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
102bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
103bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
104bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
105742fe5e2SPeter Brune   if (fas->next)      ierr = SNESReset(fas->next);CHKERRQ(ierr);
10622c1e704SPeter Brune 
107421d9b32SPeter Brune   PetscFunctionReturn(0);
108421d9b32SPeter Brune }
109421d9b32SPeter Brune 
110421d9b32SPeter Brune #undef __FUNCT__
111421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
112421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
113421d9b32SPeter Brune {
114421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
115742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
116421d9b32SPeter Brune 
117421d9b32SPeter Brune   PetscFunctionBegin;
118421d9b32SPeter Brune   /* recursively resets and then destroys */
11979d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
120421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
121421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
122ee1fd11aSPeter Brune 
123421d9b32SPeter Brune   PetscFunctionReturn(0);
124421d9b32SPeter Brune }
125421d9b32SPeter Brune 
126421d9b32SPeter Brune #undef __FUNCT__
127421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
128421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
129421d9b32SPeter Brune {
13048bfdf8aSPeter Brune   SNES_FAS                    *fas = (SNES_FAS *) snes->data;
131421d9b32SPeter Brune   PetscErrorCode              ierr;
132efe1f98aSPeter Brune   VecScatter                  injscatter;
133d1adcc6fSPeter Brune   PetscInt                    dm_levels;
1343dccd265SPeter Brune   Vec                         vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
135ab8d36c9SPeter Brune   SNES                        next;
136ab8d36c9SPeter Brune   PetscBool                   isFine;
137*f89ba88eSPeter Brune   SNESLineSearch              linesearch;
138*f89ba88eSPeter Brune   SNESLineSearch              slinesearch;
139*f89ba88eSPeter Brune   void                        *lsprectx,*lspostctx;
140*f89ba88eSPeter Brune   SNESLineSearchPreCheckFunc  lsprefunc;
141*f89ba88eSPeter Brune   SNESLineSearchPostCheckFunc lspostfunc;
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) {
173e4ed7901SPeter 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   }
182534ebe21SPeter Brune   if (!fas->smoothu && fas->level != fas->levels - 1) {
183534ebe21SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
184534ebe21SPeter Brune   }
185ab8d36c9SPeter Brune 
18679d9a41aSPeter Brune   if (snes->dm) {
187ab8d36c9SPeter Brune     /* set the smoother DMs properly */
188ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
189ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
19079d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
191ab8d36c9SPeter Brune     if (next) {
19279d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
193ab8d36c9SPeter Brune       if (!next->dm) {
194ab8d36c9SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr);
195ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
19679d9a41aSPeter Brune       }
19779d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
19879d9a41aSPeter Brune       if (!fas->interpolate) {
199ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
200bccf9bb3SJed Brown         if (!fas->restrct) {
201bccf9bb3SJed Brown           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
20279d9a41aSPeter Brune           fas->restrct = fas->interpolate;
20379d9a41aSPeter Brune         }
204bccf9bb3SJed Brown       }
20579d9a41aSPeter Brune       /* set the injection from the DM */
20679d9a41aSPeter Brune       if (!fas->inject) {
207ab8d36c9SPeter Brune         ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
20879d9a41aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
20979d9a41aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
21079d9a41aSPeter Brune       }
21179d9a41aSPeter Brune     }
21279d9a41aSPeter Brune   }
21379d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
21479d9a41aSPeter Brune   if (fas->galerkin) {
215ab8d36c9SPeter Brune     if (next)
216ab8d36c9SPeter Brune       ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
217ab8d36c9SPeter Brune     if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
218ab8d36c9SPeter Brune     if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
21979d9a41aSPeter Brune   }
22079d9a41aSPeter Brune 
221534ebe21SPeter Brune   /* sets the down (pre) smoother's default norm and sets it from options */
222534ebe21SPeter Brune   if(fas->smoothd){
223534ebe21SPeter Brune     if (fas->level == 0) {
224534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothd, SNES_NORM_NONE);CHKERRQ(ierr);
225534ebe21SPeter Brune      } else {
226534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothd, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
227534ebe21SPeter Brune     }
228534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothd);CHKERRQ(ierr);
229*f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
230*f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(fas->smoothd,&slinesearch);CHKERRQ(ierr);
231*f89ba88eSPeter Brune     ierr = SNESLineSearchGetPreCheck(linesearch,&lsprefunc,&lsprectx);CHKERRQ(ierr);
232*f89ba88eSPeter Brune     ierr = SNESLineSearchGetPostCheck(linesearch,&lspostfunc,&lspostctx);CHKERRQ(ierr);
233*f89ba88eSPeter Brune     ierr = SNESLineSearchSetPreCheck(slinesearch,lsprefunc,lsprectx);CHKERRQ(ierr);
234*f89ba88eSPeter Brune     ierr = SNESLineSearchSetPostCheck(slinesearch,lspostfunc,lspostctx);CHKERRQ(ierr);
235*f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
236534ebe21SPeter Brune   }
237534ebe21SPeter Brune 
238534ebe21SPeter Brune   /* sets the up (post) smoother's default norm and sets it from options */
239534ebe21SPeter Brune   if(fas->smoothu){
240534ebe21SPeter Brune     if (fas->level != fas->levels - 1) {
241534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothu, SNES_NORM_NONE);CHKERRQ(ierr);
242534ebe21SPeter Brune     } else {
243534ebe21SPeter Brune       ierr = SNESSetNormType(fas->smoothu, SNES_NORM_FINAL_ONLY);CHKERRQ(ierr);
244534ebe21SPeter Brune     }
245534ebe21SPeter Brune     ierr = SNESSetFromOptions(fas->smoothu);CHKERRQ(ierr);
246*f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
247*f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(fas->smoothu,&slinesearch);CHKERRQ(ierr);
248*f89ba88eSPeter Brune     ierr = SNESLineSearchGetPreCheck(linesearch,&lsprefunc,&lsprectx);CHKERRQ(ierr);
249*f89ba88eSPeter Brune     ierr = SNESLineSearchGetPostCheck(linesearch,&lspostfunc,&lspostctx);CHKERRQ(ierr);
250*f89ba88eSPeter Brune     ierr = SNESLineSearchSetPreCheck(slinesearch,lsprefunc,lsprectx);CHKERRQ(ierr);
251*f89ba88eSPeter Brune     ierr = SNESLineSearchSetPostCheck(slinesearch,lspostfunc,lspostctx);CHKERRQ(ierr);
252*f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
253534ebe21SPeter Brune   }
254d06165b7SPeter Brune 
255ab8d36c9SPeter Brune   if (next) {
25679d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
257ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
258ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
259*f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(snes,&linesearch);CHKERRQ(ierr);
260*f89ba88eSPeter Brune     ierr = SNESGetSNESLineSearch(fas->next,&slinesearch);CHKERRQ(ierr);
261*f89ba88eSPeter Brune     ierr = SNESLineSearchGetPreCheck(linesearch,&lsprefunc,&lsprectx);CHKERRQ(ierr);
262*f89ba88eSPeter Brune     ierr = SNESLineSearchGetPostCheck(linesearch,&lspostfunc,&lspostctx);CHKERRQ(ierr);
263*f89ba88eSPeter Brune     ierr = SNESLineSearchSetPreCheck(slinesearch,lsprefunc,lsprectx);CHKERRQ(ierr);
264*f89ba88eSPeter Brune     ierr = SNESLineSearchSetPostCheck(slinesearch,lspostfunc,lspostctx);CHKERRQ(ierr);
265*f89ba88eSPeter Brune     ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)linesearch, (PetscObject)slinesearch);CHKERRQ(ierr);
266ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
26779d9a41aSPeter Brune   }
2686273346dSPeter Brune   /* setup FAS work vectors */
2696273346dSPeter Brune   if (fas->galerkin) {
2706273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
2716273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
2726273346dSPeter Brune   }
273421d9b32SPeter Brune   PetscFunctionReturn(0);
274421d9b32SPeter Brune }
275421d9b32SPeter Brune 
276421d9b32SPeter Brune #undef __FUNCT__
277421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
278421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
279421d9b32SPeter Brune {
280ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
281ee78dd50SPeter Brune   PetscInt       levels = 1;
2824d26bfa5SPeter Brune   PetscBool      flg = PETSC_FALSE, upflg = PETSC_FALSE, downflg = PETSC_FALSE, monflg = PETSC_FALSE, galerkinflg = PETSC_FALSE;
283421d9b32SPeter Brune   PetscErrorCode ierr;
284ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
28507144faaSPeter Brune   SNESFASType    fastype;
286fde0ff24SPeter Brune   const char     *optionsprefix;
287f1c6b773SPeter Brune   SNESLineSearch linesearch;
28866585501SPeter Brune   PetscInt       m, n_up, n_down;
289ab8d36c9SPeter Brune   SNES           next;
290ab8d36c9SPeter Brune   PetscBool      isFine;
291421d9b32SPeter Brune 
292421d9b32SPeter Brune   PetscFunctionBegin;
293ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
294c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
295ee78dd50SPeter Brune 
296ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
297ab8d36c9SPeter Brune   if (isFine) {
298ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
299c732cbdbSBarry Smith     if (!flg && snes->dm) {
300c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
301c732cbdbSBarry Smith       levels++;
302d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
303c732cbdbSBarry Smith     }
304ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
30507144faaSPeter Brune     fastype = fas->fastype;
30607144faaSPeter Brune     ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
30707144faaSPeter Brune     if (flg) {
30807144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
30907144faaSPeter Brune     }
310ee78dd50SPeter Brune 
311fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
312ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
313ab8d36c9SPeter Brune     if (flg) {
314ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
315fde0ff24SPeter Brune     }
316fde0ff24SPeter Brune 
317ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
318ab8d36c9SPeter Brune     if (flg) {
319ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
320ab8d36c9SPeter Brune     }
321ee78dd50SPeter Brune 
32266585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",fas->max_up_it,&n_up,&upflg);CHKERRQ(ierr);
323162d76ddSPeter Brune 
32466585501SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",fas->max_down_it,&n_down,&downflg);CHKERRQ(ierr);
325162d76ddSPeter Brune 
326c8c899caSPeter Brune     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
327c8c899caSPeter Brune     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
328ab8d36c9SPeter Brune   }
329ee78dd50SPeter Brune 
330421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
3318cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
332162d76ddSPeter Brune   if (upflg) {
33366585501SPeter Brune     ierr = SNESFASSetNumberSmoothUp(snes,n_up);CHKERRQ(ierr);
334162d76ddSPeter Brune   }
335162d76ddSPeter Brune   if (downflg) {
33666585501SPeter Brune     ierr = SNESFASSetNumberSmoothDown(snes,n_down);CHKERRQ(ierr);
337162d76ddSPeter Brune   }
338eff52c0eSPeter Brune 
3399e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
3409e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
3419e764e56SPeter Brune     if (!snes->linesearch) {
342f1c6b773SPeter Brune       ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
3431a4f838cSPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNESLINESEARCHL2);CHKERRQ(ierr);
3449e764e56SPeter Brune     }
3459e764e56SPeter Brune   }
3469e764e56SPeter Brune 
347ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
348ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
349ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
350421d9b32SPeter Brune   PetscFunctionReturn(0);
351421d9b32SPeter Brune }
352421d9b32SPeter Brune 
353421d9b32SPeter Brune #undef __FUNCT__
354421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
355421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
356421d9b32SPeter Brune {
357421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
358ab8d36c9SPeter Brune   PetscBool      isFine, iascii;
359ab8d36c9SPeter Brune   PetscInt       i;
360421d9b32SPeter Brune   PetscErrorCode ierr;
361ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
362421d9b32SPeter Brune 
363421d9b32SPeter Brune   PetscFunctionBegin;
364ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
365ab8d36c9SPeter Brune   if (isFine) {
366251f4c67SDmitry Karpeev     ierr = PetscObjectTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
367421d9b32SPeter Brune     if (iascii) {
368ab8d36c9SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
369ab8d36c9SPeter Brune       if (fas->galerkin) {
370ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
371421d9b32SPeter Brune       } else {
372ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
373421d9b32SPeter Brune       }
374ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
375ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
376ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
377ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
378ab8d36c9SPeter Brune         if (!i) {
379ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
380421d9b32SPeter Brune         } else {
381ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
382421d9b32SPeter Brune         }
383ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
384ab8d36c9SPeter Brune         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
385ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
386ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
387ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
388ab8d36c9SPeter Brune         } else if (i) {
389ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
390ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
391ab8d36c9SPeter Brune           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
392ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
393ab8d36c9SPeter Brune         }
394ab8d36c9SPeter Brune       }
395421d9b32SPeter Brune     } else {
396421d9b32SPeter Brune       SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
397421d9b32SPeter Brune     }
398ab8d36c9SPeter Brune   }
399421d9b32SPeter Brune   PetscFunctionReturn(0);
400421d9b32SPeter Brune }
401421d9b32SPeter Brune 
402421d9b32SPeter Brune #undef __FUNCT__
40391f99d7cSPeter Brune #define __FUNCT__ "SNESFASDownSmooth_Private"
40439bd7f45SPeter Brune /*
40539bd7f45SPeter Brune Defines the action of the downsmoother
40639bd7f45SPeter Brune  */
40791f99d7cSPeter Brune PetscErrorCode SNESFASDownSmooth_Private(SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm)
408b9c2fdf1SPeter Brune {
40939bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
410742fe5e2SPeter Brune   SNESConvergedReason reason;
411ab8d36c9SPeter Brune   Vec                 FPC;
412ab8d36c9SPeter Brune   SNES                smoothd;
413421d9b32SPeter Brune   PetscFunctionBegin;
414ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
415e4ed7901SPeter Brune   ierr = SNESSetInitialFunction(smoothd, F);CHKERRQ(ierr);
416e4ed7901SPeter Brune   ierr = SNESSetInitialFunctionNorm(smoothd, *fnorm);CHKERRQ(ierr);
417ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
418742fe5e2SPeter Brune   /* check convergence reason for the smoother */
419ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
420e70c42e5SPeter Brune   if (reason < 0 && !(reason == SNES_DIVERGED_MAX_IT || reason == SNES_DIVERGED_LOCAL_MIN)) {
421742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
422742fe5e2SPeter Brune     PetscFunctionReturn(0);
423742fe5e2SPeter Brune   }
424ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
4254b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
426b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothd, fnorm);CHKERRQ(ierr);
42739bd7f45SPeter Brune   PetscFunctionReturn(0);
42839bd7f45SPeter Brune }
42939bd7f45SPeter Brune 
43039bd7f45SPeter Brune 
43139bd7f45SPeter Brune #undef __FUNCT__
43291f99d7cSPeter Brune #define __FUNCT__ "SNESFASUpSmooth_Private"
43339bd7f45SPeter Brune /*
43407144faaSPeter Brune Defines the action of the upsmoother
43539bd7f45SPeter Brune  */
43691f99d7cSPeter Brune PetscErrorCode SNESFASUpSmooth_Private (SNES snes, Vec B, Vec X, Vec F, PetscReal *fnorm) {
43739bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
43839bd7f45SPeter Brune   SNESConvergedReason reason;
439ab8d36c9SPeter Brune   Vec                 FPC;
440ab8d36c9SPeter Brune   SNES                smoothu;
44139bd7f45SPeter Brune   PetscFunctionBegin;
442ab8d36c9SPeter Brune 
443ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
444ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
44539bd7f45SPeter Brune   /* check convergence reason for the smoother */
446ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
44739bd7f45SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
44839bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
44939bd7f45SPeter Brune     PetscFunctionReturn(0);
45039bd7f45SPeter Brune   }
451ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
4524b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
453b9c2fdf1SPeter Brune   ierr = SNESGetFunctionNorm(smoothu, fnorm);CHKERRQ(ierr);
45439bd7f45SPeter Brune   PetscFunctionReturn(0);
45539bd7f45SPeter Brune }
45639bd7f45SPeter Brune 
45739bd7f45SPeter Brune #undef __FUNCT__
458938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
459938e4a01SJed Brown /*@
460938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
461938e4a01SJed Brown 
462938e4a01SJed Brown    Collective
463938e4a01SJed Brown 
464938e4a01SJed Brown    Input Arguments:
465938e4a01SJed Brown .  snes - SNESFAS
466938e4a01SJed Brown 
467938e4a01SJed Brown    Output Arguments:
468938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
469938e4a01SJed Brown 
470938e4a01SJed Brown    Level: developer
471938e4a01SJed Brown 
472938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
473938e4a01SJed Brown @*/
474938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
475938e4a01SJed Brown {
476938e4a01SJed Brown   PetscErrorCode ierr;
477938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
478938e4a01SJed Brown 
479938e4a01SJed Brown   PetscFunctionBegin;
480938e4a01SJed Brown   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
481938e4a01SJed Brown   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
482938e4a01SJed Brown   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
483938e4a01SJed Brown   PetscFunctionReturn(0);
484938e4a01SJed Brown }
485938e4a01SJed Brown 
486e9923e8dSJed Brown #undef __FUNCT__
487e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
488e9923e8dSJed Brown /*@
489e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
490e9923e8dSJed Brown 
491e9923e8dSJed Brown    Collective
492e9923e8dSJed Brown 
493e9923e8dSJed Brown    Input Arguments:
494e9923e8dSJed Brown +  fine - SNES from which to restrict
495e9923e8dSJed Brown -  Xfine - vector to restrict
496e9923e8dSJed Brown 
497e9923e8dSJed Brown    Output Arguments:
498e9923e8dSJed Brown .  Xcoarse - result of restriction
499e9923e8dSJed Brown 
500e9923e8dSJed Brown    Level: developer
501e9923e8dSJed Brown 
502e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
503e9923e8dSJed Brown @*/
504e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
505e9923e8dSJed Brown {
506e9923e8dSJed Brown   PetscErrorCode ierr;
507e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
508e9923e8dSJed Brown 
509e9923e8dSJed Brown   PetscFunctionBegin;
510e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
511e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
512e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
513e9923e8dSJed Brown   if (fas->inject) {
514e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
515e9923e8dSJed Brown   } else {
516e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
517e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
518e9923e8dSJed Brown   }
519e9923e8dSJed Brown   PetscFunctionReturn(0);
520e9923e8dSJed Brown }
521e9923e8dSJed Brown 
522e9923e8dSJed Brown #undef __FUNCT__
5238c40d5fbSBarry Smith #define __FUNCT__ "SNESFASCoarseCorrection"
52439bd7f45SPeter Brune /*
52539bd7f45SPeter Brune 
52639bd7f45SPeter Brune Performs the FAS coarse correction as:
52739bd7f45SPeter Brune 
52839bd7f45SPeter Brune fine problem: F(x) = 0
52939bd7f45SPeter Brune coarse problem: F^c(x) = b^c
53039bd7f45SPeter Brune 
53139bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
53239bd7f45SPeter Brune 
53339bd7f45SPeter Brune  */
5348c40d5fbSBarry Smith PetscErrorCode SNESFASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
53539bd7f45SPeter Brune   PetscErrorCode      ierr;
53639bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
53739bd7f45SPeter Brune   SNESConvergedReason reason;
538ab8d36c9SPeter Brune   SNES                next;
539ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
54039bd7f45SPeter Brune   PetscFunctionBegin;
541ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
542ab8d36c9SPeter Brune   if (next) {
543ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
544ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
545ab8d36c9SPeter Brune 
546ab8d36c9SPeter Brune     X_c  = next->vec_sol;
547ab8d36c9SPeter Brune     Xo_c = next->work[0];
548ab8d36c9SPeter Brune     F_c  = next->vec_func;
549ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
550efe1f98aSPeter Brune 
551938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
552293a7e31SPeter Brune     /* restrict the defect */
553ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
554b9c2fdf1SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - R(F(x) - b) */
555ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
556e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
557b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
558e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
559ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
560ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
561c90fad12SPeter Brune 
562c90fad12SPeter Brune     /* recurse to the next level */
563e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
564ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
565ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
566742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
567742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
568742fe5e2SPeter Brune       PetscFunctionReturn(0);
569742fe5e2SPeter Brune     }
570fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
571fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
572ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
573293a7e31SPeter Brune   }
57439bd7f45SPeter Brune   PetscFunctionReturn(0);
57539bd7f45SPeter Brune }
57639bd7f45SPeter Brune 
57739bd7f45SPeter Brune #undef __FUNCT__
5782cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Additive"
57939bd7f45SPeter Brune /*
58039bd7f45SPeter Brune 
58139bd7f45SPeter Brune The additive cycle looks like:
58239bd7f45SPeter Brune 
58307144faaSPeter Brune xhat = x
58407144faaSPeter Brune xhat = dS(x, b)
58507144faaSPeter Brune x = coarsecorrection(xhat, b_d)
58607144faaSPeter Brune x = x + nu*(xhat - x);
58739bd7f45SPeter Brune (optional) x = uS(x, b)
58839bd7f45SPeter Brune 
58939bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
59039bd7f45SPeter Brune 
59139bd7f45SPeter Brune  */
59291f99d7cSPeter Brune PetscErrorCode SNESFASCycle_Additive(SNES snes, Vec X) {
59307144faaSPeter Brune   Vec                 F, B, Xhat;
59422c1e704SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
59539bd7f45SPeter Brune   PetscErrorCode      ierr;
59607144faaSPeter Brune   SNESConvergedReason reason;
59722c1e704SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
59822c1e704SPeter Brune   PetscBool           lssuccess;
599ab8d36c9SPeter Brune   SNES                next;
600ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
60139bd7f45SPeter Brune   PetscFunctionBegin;
602ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
60339bd7f45SPeter Brune   F = snes->vec_func;
60439bd7f45SPeter Brune   B = snes->vec_rhs;
605e7f468e7SPeter Brune   Xhat = snes->work[1];
60607144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
60707144faaSPeter Brune   /* recurse first */
608ab8d36c9SPeter Brune   if (next) {
609ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
610ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
61107144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
612c2a02606SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
613ab8d36c9SPeter Brune     X_c  = next->vec_sol;
614ab8d36c9SPeter Brune     Xo_c = next->work[0];
615ab8d36c9SPeter Brune     F_c  = next->vec_func;
616ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
61739bd7f45SPeter Brune 
618938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
61907144faaSPeter Brune     /* restrict the defect */
620ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
62107144faaSPeter Brune 
62207144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
623ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
624e4ed7901SPeter Brune     ierr = VecCopy(B_c, X_c);CHKERRQ(ierr);
625b9c2fdf1SPeter Brune     ierr = VecCopy(F_c, B_c);CHKERRQ(ierr);
626e4ed7901SPeter Brune     ierr = VecCopy(X_c, F_c);CHKERRQ(ierr);
62707144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
62807144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
62907144faaSPeter Brune 
63007144faaSPeter Brune     /* recurse */
631e4ed7901SPeter Brune     ierr = SNESSetInitialFunction(next, F_c);CHKERRQ(ierr);
632ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
63307144faaSPeter Brune 
63407144faaSPeter Brune     /* smooth on this level */
63591f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &fnorm);CHKERRQ(ierr);
63607144faaSPeter Brune 
637ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
63807144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
63907144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
64007144faaSPeter Brune       PetscFunctionReturn(0);
64107144faaSPeter Brune     }
64207144faaSPeter Brune 
64307144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
644c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
645ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
64607144faaSPeter Brune 
647ddebd997SPeter Brune     /* additive correction of the coarse direction*/
648f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
649f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
6509e764e56SPeter Brune     if (!lssuccess) {
6519e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6529e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6539e764e56SPeter Brune         PetscFunctionReturn(0);
6549e764e56SPeter Brune       }
6559e764e56SPeter Brune     }
656b9c2fdf1SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &snes->norm, &ynorm);CHKERRQ(ierr);
65707144faaSPeter Brune   } else {
65891f99d7cSPeter Brune     ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
65907144faaSPeter Brune   }
66039bd7f45SPeter Brune   PetscFunctionReturn(0);
66139bd7f45SPeter Brune }
66239bd7f45SPeter Brune 
66339bd7f45SPeter Brune #undef __FUNCT__
6642cf9d1e8SPeter Brune #define __FUNCT__ "SNESFASCycle_Multiplicative"
66539bd7f45SPeter Brune /*
66639bd7f45SPeter Brune 
66739bd7f45SPeter Brune Defines the FAS cycle as:
66839bd7f45SPeter Brune 
66939bd7f45SPeter Brune fine problem: F(x) = 0
67039bd7f45SPeter Brune coarse problem: F^c(x) = b^c
67139bd7f45SPeter Brune 
67239bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
67339bd7f45SPeter Brune 
67439bd7f45SPeter Brune correction:
67539bd7f45SPeter Brune 
67639bd7f45SPeter Brune x = x + I(x^c - Rx)
67739bd7f45SPeter Brune 
67839bd7f45SPeter Brune  */
67991f99d7cSPeter Brune PetscErrorCode SNESFASCycle_Multiplicative(SNES snes, Vec X) {
68039bd7f45SPeter Brune 
68139bd7f45SPeter Brune   PetscErrorCode      ierr;
68239bd7f45SPeter Brune   Vec                 F,B;
68339bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
68439bd7f45SPeter Brune 
68539bd7f45SPeter Brune   PetscFunctionBegin;
68639bd7f45SPeter Brune   F = snes->vec_func;
68739bd7f45SPeter Brune   B = snes->vec_rhs;
68839bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
68991f99d7cSPeter Brune   ierr = SNESFASDownSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
690c90fad12SPeter Brune   if (fas->level != 0) {
6918c40d5fbSBarry Smith     ierr = SNESFASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
69291f99d7cSPeter Brune     ierr = SNESFASUpSmooth_Private(snes, B, X, F, &snes->norm);CHKERRQ(ierr);
693fe6f9142SPeter Brune   }
694fa9694d7SPeter Brune   PetscFunctionReturn(0);
695421d9b32SPeter Brune }
696421d9b32SPeter Brune 
697421d9b32SPeter Brune #undef __FUNCT__
698421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
699421d9b32SPeter Brune 
700421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
701421d9b32SPeter Brune {
702fa9694d7SPeter Brune   PetscErrorCode ierr;
703fe6f9142SPeter Brune   PetscInt       i, maxits;
704ddb5aff1SPeter Brune   Vec            X, F;
705fe6f9142SPeter Brune   PetscReal      fnorm;
706b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
707b17ce1afSJed Brown   DM             dm;
708e70c42e5SPeter Brune   PetscBool      isFine;
709b17ce1afSJed Brown 
710421d9b32SPeter Brune   PetscFunctionBegin;
711fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
712fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
713fa9694d7SPeter Brune   X = snes->vec_sol;
714f5a6d4f9SBarry Smith   F = snes->vec_func;
715293a7e31SPeter Brune 
71618a66777SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
717293a7e31SPeter Brune   /*norm setup */
718fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
719fe6f9142SPeter Brune   snes->iter = 0;
720fe6f9142SPeter Brune   snes->norm = 0.;
721fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
722e4ed7901SPeter Brune   if (!snes->vec_func_init_set) {
723fe6f9142SPeter Brune     ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
724fe6f9142SPeter Brune     if (snes->domainerror) {
725fe6f9142SPeter Brune       snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
726fe6f9142SPeter Brune       PetscFunctionReturn(0);
727fe6f9142SPeter Brune     }
728e4ed7901SPeter Brune   } else {
729e4ed7901SPeter Brune     snes->vec_func_init_set = PETSC_FALSE;
730e4ed7901SPeter Brune   }
731e4ed7901SPeter Brune 
732e4ed7901SPeter Brune   if (!snes->norm_init_set) {
733fe6f9142SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
734fe6f9142SPeter Brune     if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
735fe6f9142SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
736e4ed7901SPeter Brune   } else {
737e4ed7901SPeter Brune     fnorm = snes->norm_init;
738e4ed7901SPeter Brune     snes->norm_init_set = PETSC_FALSE;
739e4ed7901SPeter Brune   }
740e4ed7901SPeter Brune 
741fe6f9142SPeter Brune   snes->norm = fnorm;
742fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
743fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
744fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
745fe6f9142SPeter Brune 
746fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
747fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
748fe6f9142SPeter Brune   /* test convergence */
749fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
750fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
751e4ed7901SPeter Brune 
752b17ce1afSJed Brown 
753b9c2fdf1SPeter Brune   if (isFine) {
754b9c2fdf1SPeter Brune     /* propagate scale-dependent data up the hierarchy */
755b17ce1afSJed Brown     ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
756b17ce1afSJed Brown     for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
757b17ce1afSJed Brown       DM dmcoarse;
758b17ce1afSJed Brown       ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
759b17ce1afSJed Brown       ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
760b17ce1afSJed Brown       dm = dmcoarse;
761b17ce1afSJed Brown     }
762b9c2fdf1SPeter Brune   }
763b17ce1afSJed Brown 
764fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
765fe6f9142SPeter Brune     /* Call general purpose update function */
766646217ecSPeter Brune 
767fe6f9142SPeter Brune     if (snes->ops->update) {
768fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
769fe6f9142SPeter Brune     }
77007144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
77191f99d7cSPeter Brune       ierr = SNESFASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
77207144faaSPeter Brune     } else {
77391f99d7cSPeter Brune       ierr = SNESFASCycle_Additive(snes, X);CHKERRQ(ierr);
77407144faaSPeter Brune     }
775742fe5e2SPeter Brune 
776742fe5e2SPeter Brune     /* check for FAS cycle divergence */
777742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
778742fe5e2SPeter Brune       PetscFunctionReturn(0);
779742fe5e2SPeter Brune     }
780b9c2fdf1SPeter Brune 
781c90fad12SPeter Brune     /* Monitor convergence */
782c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
783c90fad12SPeter Brune     snes->iter = i+1;
784c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
785c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
786c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
787c90fad12SPeter Brune     /* Test for convergence */
78866585501SPeter Brune     if (isFine) {
789b9c2fdf1SPeter Brune       ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,snes->norm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
790c90fad12SPeter Brune       if (snes->reason) break;
791fe6f9142SPeter Brune     }
79266585501SPeter Brune   }
793fe6f9142SPeter Brune   if (i == maxits) {
794fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
795fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
796fe6f9142SPeter Brune   }
797421d9b32SPeter Brune   PetscFunctionReturn(0);
798421d9b32SPeter Brune }
799