xref: /petsc/src/snes/impls/fas/fas.c (revision ab8d36c96b69fea5736263697ddaa972e87bcd95)
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 *);
13*ab8d36c9SPeter 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 
191fbfccc6SJed Brown The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
201fbfccc6SJed Brown 
211fbfccc6SJed Brown Level: advanced
221fbfccc6SJed Brown 
231fbfccc6SJed Brown .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
241fbfccc6SJed Brown M*/
25421d9b32SPeter Brune 
26421d9b32SPeter Brune #undef __FUNCT__
27421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
281fbfccc6SJed Brown PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes)
29421d9b32SPeter Brune {
30421d9b32SPeter Brune   SNES_FAS * fas;
31421d9b32SPeter Brune   PetscErrorCode ierr;
32421d9b32SPeter Brune 
33421d9b32SPeter Brune   PetscFunctionBegin;
34421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
35421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
36421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
37421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
38421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
39421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
40421d9b32SPeter Brune 
41ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
42ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
43ed020824SBarry Smith 
4488976e71SPeter Brune   if (!snes->tolerancesset) {
450e444f03SPeter Brune     snes->max_funcs = 30000;
460e444f03SPeter Brune     snes->max_its   = 10000;
4788976e71SPeter Brune   }
480e444f03SPeter Brune 
49421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
50421d9b32SPeter Brune   snes->data                  = (void*) fas;
51421d9b32SPeter Brune   fas->level                  = 0;
52293a7e31SPeter Brune   fas->levels                 = 1;
53ee78dd50SPeter Brune   fas->n_cycles               = 1;
54ee78dd50SPeter Brune   fas->max_up_it              = 1;
55ee78dd50SPeter Brune   fas->max_down_it            = 1;
56*ab8d36c9SPeter Brune   fas->smoothu                = PETSC_NULL;
57*ab8d36c9SPeter Brune   fas->smoothd                = PETSC_NULL;
58421d9b32SPeter Brune   fas->next                   = PETSC_NULL;
596273346dSPeter Brune   fas->previous               = PETSC_NULL;
60*ab8d36c9SPeter Brune   fas->fine                   = snes;
61421d9b32SPeter Brune   fas->interpolate            = PETSC_NULL;
62421d9b32SPeter Brune   fas->restrct                = PETSC_NULL;
63efe1f98aSPeter Brune   fas->inject                 = PETSC_NULL;
64cc05f883SPeter Brune   fas->monitor                = PETSC_NULL;
65cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
66ddebd997SPeter Brune   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
67efe1f98aSPeter Brune   PetscFunctionReturn(0);
68efe1f98aSPeter Brune }
69efe1f98aSPeter Brune 
70421d9b32SPeter Brune #undef __FUNCT__
71421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
72421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
73421d9b32SPeter Brune {
7477df8cc4SPeter Brune   PetscErrorCode ierr = 0;
75421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
76421d9b32SPeter Brune 
77421d9b32SPeter Brune   PetscFunctionBegin;
78*ab8d36c9SPeter Brune   if (fas->smoothu != fas->smoothd) {
79*ab8d36c9SPeter Brune     ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
80*ab8d36c9SPeter Brune     ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
81*ab8d36c9SPeter Brune   } else {
82*ab8d36c9SPeter Brune     ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
83*ab8d36c9SPeter Brune     fas->smoothu = PETSC_NULL;
84*ab8d36c9SPeter Brune   }
853dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
86bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
87bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
88bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
89742fe5e2SPeter Brune   if (fas->next)      ierr = SNESReset(fas->next);CHKERRQ(ierr);
9022c1e704SPeter Brune 
91421d9b32SPeter Brune   PetscFunctionReturn(0);
92421d9b32SPeter Brune }
93421d9b32SPeter Brune 
94421d9b32SPeter Brune #undef __FUNCT__
95421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
96421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
97421d9b32SPeter Brune {
98421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
99742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
100421d9b32SPeter Brune 
101421d9b32SPeter Brune   PetscFunctionBegin;
102421d9b32SPeter Brune   /* recursively resets and then destroys */
10379d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
104421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
105421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
106ee1fd11aSPeter Brune 
107421d9b32SPeter Brune   PetscFunctionReturn(0);
108421d9b32SPeter Brune }
109421d9b32SPeter Brune 
110421d9b32SPeter Brune #undef __FUNCT__
111421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
112421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
113421d9b32SPeter Brune {
11448bfdf8aSPeter Brune   SNES_FAS                *fas = (SNES_FAS *) snes->data;
115421d9b32SPeter Brune   PetscErrorCode          ierr;
116efe1f98aSPeter Brune   VecScatter              injscatter;
117d1adcc6fSPeter Brune   PetscInt                dm_levels;
1183dccd265SPeter Brune   Vec                     vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
119*ab8d36c9SPeter Brune   SNES                    next;
120*ab8d36c9SPeter Brune   PetscBool               isFine;
121421d9b32SPeter Brune   PetscFunctionBegin;
122eff52c0eSPeter Brune 
123*ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
124*ab8d36c9SPeter Brune   if (fas->usedmfornumberoflevels && isFine) {
125d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
126d1adcc6fSPeter Brune     dm_levels++;
127cc05f883SPeter Brune     if (dm_levels > fas->levels) {
1282e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
1293dccd265SPeter Brune       vec_sol = snes->vec_sol;
1303dccd265SPeter Brune       vec_func = snes->vec_func;
1313dccd265SPeter Brune       vec_sol_update = snes->vec_sol_update;
1323dccd265SPeter Brune       vec_rhs = snes->vec_rhs;
1333dccd265SPeter Brune       snes->vec_sol = PETSC_NULL;
1343dccd265SPeter Brune       snes->vec_func = PETSC_NULL;
1353dccd265SPeter Brune       snes->vec_sol_update = PETSC_NULL;
1363dccd265SPeter Brune       snes->vec_rhs = PETSC_NULL;
1373dccd265SPeter Brune 
1383dccd265SPeter Brune       /* reset the number of levels */
139d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
140cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
1413dccd265SPeter Brune 
1423dccd265SPeter Brune       snes->vec_sol = vec_sol;
1433dccd265SPeter Brune       snes->vec_func = vec_func;
1443dccd265SPeter Brune       snes->vec_rhs = vec_rhs;
1453dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
146d1adcc6fSPeter Brune     }
147d1adcc6fSPeter Brune   }
148*ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
149*ab8d36c9SPeter Brune   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
1503dccd265SPeter Brune 
15107144faaSPeter Brune   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
152fde0ff24SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
15307144faaSPeter Brune   } else {
154ddebd997SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
15507144faaSPeter Brune   }
156cc05f883SPeter Brune 
157*ab8d36c9SPeter Brune   /* set up the smoothers if they haven't already been set up */
158*ab8d36c9SPeter Brune   if (!fas->smoothd) {
159*ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
160*ab8d36c9SPeter Brune   }
161*ab8d36c9SPeter Brune 
16279d9a41aSPeter Brune   if (snes->dm) {
163*ab8d36c9SPeter Brune     /* set the smoother DMs properly */
164*ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
165*ab8d36c9SPeter Brune     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
16679d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
167*ab8d36c9SPeter Brune     if (next) {
16879d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
169*ab8d36c9SPeter Brune       if (!next->dm) {
170*ab8d36c9SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr);
171*ab8d36c9SPeter Brune         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
17279d9a41aSPeter Brune       }
17379d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
17479d9a41aSPeter Brune       if (!fas->interpolate) {
175*ab8d36c9SPeter Brune         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
176bccf9bb3SJed Brown         if (!fas->restrct) {
177bccf9bb3SJed Brown           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
17879d9a41aSPeter Brune           fas->restrct = fas->interpolate;
17979d9a41aSPeter Brune         }
180bccf9bb3SJed Brown       }
18179d9a41aSPeter Brune       /* set the injection from the DM */
18279d9a41aSPeter Brune       if (!fas->inject) {
183*ab8d36c9SPeter Brune         ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
18479d9a41aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
18579d9a41aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
18679d9a41aSPeter Brune       }
18779d9a41aSPeter Brune     }
18879d9a41aSPeter Brune   }
18979d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
19079d9a41aSPeter Brune 
19179d9a41aSPeter Brune   if (fas->galerkin) {
192*ab8d36c9SPeter Brune     if (next)
193*ab8d36c9SPeter Brune       ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
194*ab8d36c9SPeter Brune     if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
195*ab8d36c9SPeter Brune     if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
19679d9a41aSPeter Brune   }
19779d9a41aSPeter Brune 
198*ab8d36c9SPeter Brune   if (next) {
19979d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
200*ab8d36c9SPeter Brune     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
201*ab8d36c9SPeter Brune     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
202*ab8d36c9SPeter Brune     ierr = SNESSetUp(next);CHKERRQ(ierr);
20379d9a41aSPeter Brune   }
2046273346dSPeter Brune   /* setup FAS work vectors */
2056273346dSPeter Brune   if (fas->galerkin) {
2066273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
2076273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
2086273346dSPeter Brune   }
209421d9b32SPeter Brune   PetscFunctionReturn(0);
210421d9b32SPeter Brune }
211421d9b32SPeter Brune 
212421d9b32SPeter Brune #undef __FUNCT__
213421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
214421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
215421d9b32SPeter Brune {
216ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
217ee78dd50SPeter Brune   PetscInt       levels = 1;
218*ab8d36c9SPeter Brune   PetscBool      flg, monflg, galerkinflg;
219421d9b32SPeter Brune   PetscErrorCode ierr;
220ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
22107144faaSPeter Brune   SNESFASType    fastype;
222fde0ff24SPeter Brune   const char     *optionsprefix;
223f1c6b773SPeter Brune   SNESLineSearch linesearch;
224*ab8d36c9SPeter Brune   PetscInt       m;
225*ab8d36c9SPeter Brune   SNES           next;
226*ab8d36c9SPeter Brune   PetscBool      isFine;
227421d9b32SPeter Brune 
228421d9b32SPeter Brune   PetscFunctionBegin;
229*ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
230c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
231ee78dd50SPeter Brune 
232*ab8d36c9SPeter Brune   /* number of levels -- only process most options on the finest level */
233*ab8d36c9SPeter Brune   if (isFine) {
234ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
235c732cbdbSBarry Smith     if (!flg && snes->dm) {
236c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
237c732cbdbSBarry Smith       levels++;
238d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
239c732cbdbSBarry Smith     }
240ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
24107144faaSPeter Brune     fastype = fas->fastype;
24207144faaSPeter Brune     ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
24307144faaSPeter Brune     if (flg) {
24407144faaSPeter Brune       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
24507144faaSPeter Brune     }
246ee78dd50SPeter Brune 
247fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
248*ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
249*ab8d36c9SPeter Brune     if (flg) {
250*ab8d36c9SPeter Brune       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
251fde0ff24SPeter Brune     }
252fde0ff24SPeter Brune 
253*ab8d36c9SPeter Brune     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
254*ab8d36c9SPeter Brune     if (flg) {
255*ab8d36c9SPeter Brune       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
256*ab8d36c9SPeter Brune     }
257ee78dd50SPeter Brune 
258*ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
259*ab8d36c9SPeter Brune     if (flg) {
260*ab8d36c9SPeter Brune       ierr = SNESFASSetNumberSmoothUp(snes,m);CHKERRQ(ierr);
261*ab8d36c9SPeter Brune     }
262*ab8d36c9SPeter Brune     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
263*ab8d36c9SPeter Brune     if (flg) {
264*ab8d36c9SPeter Brune       ierr = SNESFASSetNumberSmoothDown(snes,m);CHKERRQ(ierr);
265*ab8d36c9SPeter Brune     }
266646217ecSPeter Brune     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
267*ab8d36c9SPeter Brune   }
268ee78dd50SPeter Brune 
269421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
2708cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
271eff52c0eSPeter Brune 
272ee78dd50SPeter Brune   if (monflg) {
273646217ecSPeter Brune     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
274794bee33SPeter Brune     /* set the monitors for the upsmoother and downsmoother */
2752f7ea302SPeter Brune     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
276742fe5e2SPeter Brune     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
277*ab8d36c9SPeter Brune     if (fas->smoothu) ierr = SNESMonitorSet(fas->smoothu,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
278*ab8d36c9SPeter Brune     if (fas->smoothd) ierr = SNESMonitorSet(fas->smoothd,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
279d28543b3SPeter Brune   } else {
280d28543b3SPeter Brune     /* unset the monitors on the coarse levels */
281*ab8d36c9SPeter Brune     if (!isFine) {
282d28543b3SPeter Brune       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
283d28543b3SPeter Brune     }
284ee78dd50SPeter Brune   }
285ee78dd50SPeter Brune 
2869e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
2879e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
2889e764e56SPeter Brune     if (!snes->linesearch) {
289f1c6b773SPeter Brune       ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
290f1c6b773SPeter Brune       ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr);
2919e764e56SPeter Brune     }
2929e764e56SPeter Brune   }
2939e764e56SPeter Brune 
294*ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
295ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
296*ab8d36c9SPeter Brune   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
297421d9b32SPeter Brune   PetscFunctionReturn(0);
298421d9b32SPeter Brune }
299421d9b32SPeter Brune 
300421d9b32SPeter Brune #undef __FUNCT__
301421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
302421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
303421d9b32SPeter Brune {
304421d9b32SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
305*ab8d36c9SPeter Brune   PetscBool      isFine, iascii;
306*ab8d36c9SPeter Brune   PetscInt       i;
307421d9b32SPeter Brune   PetscErrorCode ierr;
308*ab8d36c9SPeter Brune   SNES           smoothu, smoothd, levelsnes;
309421d9b32SPeter Brune 
310421d9b32SPeter Brune   PetscFunctionBegin;
311*ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
312*ab8d36c9SPeter Brune   if (isFine) {
313421d9b32SPeter Brune     ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
314421d9b32SPeter Brune     if (iascii) {
315*ab8d36c9SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
316*ab8d36c9SPeter Brune       if (fas->galerkin) {
317*ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
318421d9b32SPeter Brune       } else {
319*ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
320421d9b32SPeter Brune       }
321*ab8d36c9SPeter Brune       for (i=0; i<fas->levels; i++) {
322*ab8d36c9SPeter Brune         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
323*ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
324*ab8d36c9SPeter Brune         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
325*ab8d36c9SPeter Brune         if (!i) {
326*ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
327421d9b32SPeter Brune         } else {
328*ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
329421d9b32SPeter Brune         }
330*ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
331*ab8d36c9SPeter Brune         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
332*ab8d36c9SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
333*ab8d36c9SPeter Brune         if (i && (smoothd == smoothu)) {
334*ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
335*ab8d36c9SPeter Brune         } else if (i) {
336*ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
337*ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
338*ab8d36c9SPeter Brune           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
339*ab8d36c9SPeter Brune           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
340*ab8d36c9SPeter Brune         }
341*ab8d36c9SPeter Brune       }
342421d9b32SPeter Brune     } else {
343421d9b32SPeter Brune       SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
344421d9b32SPeter Brune     }
345*ab8d36c9SPeter Brune   }
346421d9b32SPeter Brune   PetscFunctionReturn(0);
347421d9b32SPeter Brune }
348421d9b32SPeter Brune 
349421d9b32SPeter Brune #undef __FUNCT__
35039bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth"
35139bd7f45SPeter Brune /*
35239bd7f45SPeter Brune Defines the action of the downsmoother
35339bd7f45SPeter Brune  */
35439bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
35539bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
356742fe5e2SPeter Brune   SNESConvergedReason reason;
357*ab8d36c9SPeter Brune   Vec                 FPC;
358*ab8d36c9SPeter Brune   SNES                smoothd;
359421d9b32SPeter Brune   PetscFunctionBegin;
360*ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
361*ab8d36c9SPeter Brune   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
362742fe5e2SPeter Brune   /* check convergence reason for the smoother */
363*ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
364742fe5e2SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
365742fe5e2SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
366742fe5e2SPeter Brune     PetscFunctionReturn(0);
367742fe5e2SPeter Brune   }
368*ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
3694b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
37039bd7f45SPeter Brune   PetscFunctionReturn(0);
37139bd7f45SPeter Brune }
37239bd7f45SPeter Brune 
37339bd7f45SPeter Brune 
37439bd7f45SPeter Brune #undef __FUNCT__
37539bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth"
37639bd7f45SPeter Brune /*
37707144faaSPeter Brune Defines the action of the upsmoother
37839bd7f45SPeter Brune  */
37939bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
38039bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
38139bd7f45SPeter Brune   SNESConvergedReason reason;
382*ab8d36c9SPeter Brune   Vec                 FPC;
383*ab8d36c9SPeter Brune   SNES                smoothu;
38439bd7f45SPeter Brune   PetscFunctionBegin;
385*ab8d36c9SPeter Brune 
386*ab8d36c9SPeter Brune   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
387*ab8d36c9SPeter Brune   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
38839bd7f45SPeter Brune   /* check convergence reason for the smoother */
389*ab8d36c9SPeter Brune   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
39039bd7f45SPeter Brune   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
39139bd7f45SPeter Brune     snes->reason = SNES_DIVERGED_INNER;
39239bd7f45SPeter Brune     PetscFunctionReturn(0);
39339bd7f45SPeter Brune   }
394*ab8d36c9SPeter Brune   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
3954b32a720SPeter Brune   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
39639bd7f45SPeter Brune   PetscFunctionReturn(0);
39739bd7f45SPeter Brune }
39839bd7f45SPeter Brune 
39939bd7f45SPeter Brune #undef __FUNCT__
400938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
401938e4a01SJed Brown /*@
402938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
403938e4a01SJed Brown 
404938e4a01SJed Brown    Collective
405938e4a01SJed Brown 
406938e4a01SJed Brown    Input Arguments:
407938e4a01SJed Brown .  snes - SNESFAS
408938e4a01SJed Brown 
409938e4a01SJed Brown    Output Arguments:
410938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
411938e4a01SJed Brown 
412938e4a01SJed Brown    Level: developer
413938e4a01SJed Brown 
414938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
415938e4a01SJed Brown @*/
416938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
417938e4a01SJed Brown {
418938e4a01SJed Brown   PetscErrorCode ierr;
419938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
420938e4a01SJed Brown 
421938e4a01SJed Brown   PetscFunctionBegin;
422938e4a01SJed Brown   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
423938e4a01SJed Brown   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
424938e4a01SJed Brown   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
425938e4a01SJed Brown   PetscFunctionReturn(0);
426938e4a01SJed Brown }
427938e4a01SJed Brown 
428e9923e8dSJed Brown #undef __FUNCT__
429e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
430e9923e8dSJed Brown /*@
431e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
432e9923e8dSJed Brown 
433e9923e8dSJed Brown    Collective
434e9923e8dSJed Brown 
435e9923e8dSJed Brown    Input Arguments:
436e9923e8dSJed Brown +  fine - SNES from which to restrict
437e9923e8dSJed Brown -  Xfine - vector to restrict
438e9923e8dSJed Brown 
439e9923e8dSJed Brown    Output Arguments:
440e9923e8dSJed Brown .  Xcoarse - result of restriction
441e9923e8dSJed Brown 
442e9923e8dSJed Brown    Level: developer
443e9923e8dSJed Brown 
444e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
445e9923e8dSJed Brown @*/
446e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
447e9923e8dSJed Brown {
448e9923e8dSJed Brown   PetscErrorCode ierr;
449e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
450e9923e8dSJed Brown 
451e9923e8dSJed Brown   PetscFunctionBegin;
452e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
453e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
454e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
455e9923e8dSJed Brown   if (fas->inject) {
456e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
457e9923e8dSJed Brown   } else {
458e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
459e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
460e9923e8dSJed Brown   }
461e9923e8dSJed Brown   PetscFunctionReturn(0);
462e9923e8dSJed Brown }
463e9923e8dSJed Brown 
464e9923e8dSJed Brown #undef __FUNCT__
46539bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection"
46639bd7f45SPeter Brune /*
46739bd7f45SPeter Brune 
46839bd7f45SPeter Brune Performs the FAS coarse correction as:
46939bd7f45SPeter Brune 
47039bd7f45SPeter Brune fine problem: F(x) = 0
47139bd7f45SPeter Brune coarse problem: F^c(x) = b^c
47239bd7f45SPeter Brune 
47339bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
47439bd7f45SPeter Brune 
47539bd7f45SPeter Brune  */
47639a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
47739bd7f45SPeter Brune   PetscErrorCode      ierr;
47839bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
47939bd7f45SPeter Brune   SNESConvergedReason reason;
480*ab8d36c9SPeter Brune   SNES                next;
481*ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
48239bd7f45SPeter Brune   PetscFunctionBegin;
483*ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
484*ab8d36c9SPeter Brune   if (next) {
485*ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
486*ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
487*ab8d36c9SPeter Brune 
488*ab8d36c9SPeter Brune     X_c  = next->vec_sol;
489*ab8d36c9SPeter Brune     Xo_c = next->work[0];
490*ab8d36c9SPeter Brune     F_c  = next->vec_func;
491*ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
492efe1f98aSPeter Brune 
493938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
494293a7e31SPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
495293a7e31SPeter Brune 
496293a7e31SPeter Brune     /* restrict the defect */
497*ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
498293a7e31SPeter Brune 
499c90fad12SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
500*ab8d36c9SPeter Brune     next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
501*ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
502742fe5e2SPeter Brune 
503293a7e31SPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
504c90fad12SPeter Brune 
505ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
506ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
507c90fad12SPeter Brune 
508c90fad12SPeter Brune     /* recurse to the next level */
509*ab8d36c9SPeter Brune     next->vec_rhs = B_c;
510*ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
511*ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
512742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
513742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
514742fe5e2SPeter Brune       PetscFunctionReturn(0);
515742fe5e2SPeter Brune     }
516ee78dd50SPeter Brune 
517fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
518fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
519*ab8d36c9SPeter Brune     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
520293a7e31SPeter Brune   }
52139bd7f45SPeter Brune   PetscFunctionReturn(0);
52239bd7f45SPeter Brune }
52339bd7f45SPeter Brune 
52439bd7f45SPeter Brune #undef __FUNCT__
52539bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive"
52639bd7f45SPeter Brune /*
52739bd7f45SPeter Brune 
52839bd7f45SPeter Brune The additive cycle looks like:
52939bd7f45SPeter Brune 
53007144faaSPeter Brune xhat = x
53107144faaSPeter Brune xhat = dS(x, b)
53207144faaSPeter Brune x = coarsecorrection(xhat, b_d)
53307144faaSPeter Brune x = x + nu*(xhat - x);
53439bd7f45SPeter Brune (optional) x = uS(x, b)
53539bd7f45SPeter Brune 
53639bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
53739bd7f45SPeter Brune 
53839bd7f45SPeter Brune  */
53939bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
54007144faaSPeter Brune   Vec                 F, B, Xhat;
54122c1e704SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
54239bd7f45SPeter Brune   PetscErrorCode      ierr;
54307144faaSPeter Brune   SNESConvergedReason reason;
54422c1e704SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
54522c1e704SPeter Brune   PetscBool           lssuccess;
546*ab8d36c9SPeter Brune   SNES                next;
547*ab8d36c9SPeter Brune   Mat                 restrct, interpolate;
54839bd7f45SPeter Brune   PetscFunctionBegin;
549*ab8d36c9SPeter Brune   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
55039bd7f45SPeter Brune   F = snes->vec_func;
55139bd7f45SPeter Brune   B = snes->vec_rhs;
552fde0ff24SPeter Brune   Xhat = snes->work[3];
55307144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
55407144faaSPeter Brune   /* recurse first */
555*ab8d36c9SPeter Brune   if (next) {
556*ab8d36c9SPeter Brune     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
557*ab8d36c9SPeter Brune     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
55807144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
55939bd7f45SPeter Brune 
560*ab8d36c9SPeter Brune     X_c  = next->vec_sol;
561*ab8d36c9SPeter Brune     Xo_c = next->work[0];
562*ab8d36c9SPeter Brune     F_c  = next->vec_func;
563*ab8d36c9SPeter Brune     B_c  = next->vec_rhs;
56439bd7f45SPeter Brune 
565938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
56607144faaSPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
56707144faaSPeter Brune 
56807144faaSPeter Brune     /* restrict the defect */
569*ab8d36c9SPeter Brune     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
57007144faaSPeter Brune 
57107144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
572*ab8d36c9SPeter Brune     next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
573*ab8d36c9SPeter Brune     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
57407144faaSPeter Brune 
57507144faaSPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
57607144faaSPeter Brune 
57707144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
57807144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
57907144faaSPeter Brune 
58007144faaSPeter Brune     /* recurse */
581*ab8d36c9SPeter Brune     next->vec_rhs = B_c;
582*ab8d36c9SPeter Brune     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
58307144faaSPeter Brune 
58407144faaSPeter Brune     /* smooth on this level */
58507144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
58607144faaSPeter Brune 
587*ab8d36c9SPeter Brune     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
58807144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
58907144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
59007144faaSPeter Brune       PetscFunctionReturn(0);
59107144faaSPeter Brune     }
59207144faaSPeter Brune 
59307144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
594c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
595*ab8d36c9SPeter Brune     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
59607144faaSPeter Brune 
597ddebd997SPeter Brune     /* additive correction of the coarse direction*/
598ddebd997SPeter Brune     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
599ddebd997SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
600f1c6b773SPeter Brune     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
601f1c6b773SPeter Brune     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
6029e764e56SPeter Brune     if (!lssuccess) {
6039e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
6049e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
6059e764e56SPeter Brune         PetscFunctionReturn(0);
6069e764e56SPeter Brune       }
6079e764e56SPeter Brune     }
608f1c6b773SPeter Brune     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
60907144faaSPeter Brune   } else {
61007144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
61107144faaSPeter Brune   }
61239bd7f45SPeter Brune   PetscFunctionReturn(0);
61339bd7f45SPeter Brune }
61439bd7f45SPeter Brune 
61539bd7f45SPeter Brune #undef __FUNCT__
61639bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative"
61739bd7f45SPeter Brune /*
61839bd7f45SPeter Brune 
61939bd7f45SPeter Brune Defines the FAS cycle as:
62039bd7f45SPeter Brune 
62139bd7f45SPeter Brune fine problem: F(x) = 0
62239bd7f45SPeter Brune coarse problem: F^c(x) = b^c
62339bd7f45SPeter Brune 
62439bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
62539bd7f45SPeter Brune 
62639bd7f45SPeter Brune correction:
62739bd7f45SPeter Brune 
62839bd7f45SPeter Brune x = x + I(x^c - Rx)
62939bd7f45SPeter Brune 
63039bd7f45SPeter Brune  */
63139bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
63239bd7f45SPeter Brune 
63339bd7f45SPeter Brune   PetscErrorCode      ierr;
63439bd7f45SPeter Brune   Vec                 F,B;
63539bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
63639bd7f45SPeter Brune 
63739bd7f45SPeter Brune   PetscFunctionBegin;
63839bd7f45SPeter Brune   F = snes->vec_func;
63939bd7f45SPeter Brune   B = snes->vec_rhs;
64039bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
64139bd7f45SPeter Brune   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
64239bd7f45SPeter Brune 
64339a4b097SPeter Brune   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
64439bd7f45SPeter Brune 
645c90fad12SPeter Brune   if (fas->level != 0) {
64639bd7f45SPeter Brune     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
647fe6f9142SPeter Brune   }
648fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
649fa9694d7SPeter 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;
664b17ce1afSJed Brown 
665421d9b32SPeter Brune   PetscFunctionBegin;
666fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
667fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
668fa9694d7SPeter Brune   X = snes->vec_sol;
669f5a6d4f9SBarry Smith   F = snes->vec_func;
670293a7e31SPeter Brune 
671293a7e31SPeter Brune   /*norm setup */
672fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
673fe6f9142SPeter Brune   snes->iter = 0;
674fe6f9142SPeter Brune   snes->norm = 0.;
675fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
676fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
677fe6f9142SPeter Brune   if (snes->domainerror) {
678fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
679fe6f9142SPeter Brune     PetscFunctionReturn(0);
680fe6f9142SPeter Brune   }
681fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
682fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
683fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
684fe6f9142SPeter Brune   snes->norm = fnorm;
685fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
686fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
687fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
688fe6f9142SPeter Brune 
689fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
690fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
691fe6f9142SPeter Brune   /* test convergence */
692fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
693fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
694b17ce1afSJed Brown 
695b17ce1afSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
696b17ce1afSJed Brown   for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
697b17ce1afSJed Brown     DM dmcoarse;
698b17ce1afSJed Brown     ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
699b17ce1afSJed Brown     ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
700b17ce1afSJed Brown     dm = dmcoarse;
701b17ce1afSJed Brown   }
702b17ce1afSJed Brown 
703fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
704fe6f9142SPeter Brune     /* Call general purpose update function */
705646217ecSPeter Brune 
706fe6f9142SPeter Brune     if (snes->ops->update) {
707fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
708fe6f9142SPeter Brune     }
70907144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
71039bd7f45SPeter Brune       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
71107144faaSPeter Brune     } else {
71207144faaSPeter Brune       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
71307144faaSPeter Brune     }
714742fe5e2SPeter Brune 
715742fe5e2SPeter Brune     /* check for FAS cycle divergence */
716742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
717742fe5e2SPeter Brune       PetscFunctionReturn(0);
718742fe5e2SPeter Brune     }
719c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
720c90fad12SPeter Brune     /* Monitor convergence */
721c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
722c90fad12SPeter Brune     snes->iter = i+1;
723c90fad12SPeter Brune     snes->norm = fnorm;
724c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
725c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
726c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
727c90fad12SPeter Brune     /* Test for convergence */
728c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
729c90fad12SPeter Brune     if (snes->reason) break;
730fe6f9142SPeter Brune   }
731fe6f9142SPeter Brune   if (i == maxits) {
732fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
733fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
734fe6f9142SPeter Brune   }
735421d9b32SPeter Brune   PetscFunctionReturn(0);
736421d9b32SPeter Brune }
737