xref: /petsc/src/snes/impls/fas/fas.c (revision c732cbdb66a8ef1bb131466e08faf4a2cdd6079c)
1421d9b32SPeter Brune /* Defines the basic SNES object */
2421d9b32SPeter Brune #include <../src/snes/impls/fas/fasimpls.h>
3421d9b32SPeter Brune 
4421d9b32SPeter Brune /*MC
5421d9b32SPeter Brune Full Approximation Scheme nonlinear multigrid solver.
6421d9b32SPeter Brune 
7421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
8421d9b32SPeter Brune 
9421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
10421d9b32SPeter Brune M*/
11421d9b32SPeter Brune 
12421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
13421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
14421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
15421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
16421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
17421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
18421d9b32SPeter Brune 
19421d9b32SPeter Brune EXTERN_C_BEGIN
20421d9b32SPeter Brune 
21421d9b32SPeter Brune #undef __FUNCT__
22421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
23421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes)
24421d9b32SPeter Brune {
25421d9b32SPeter Brune   SNES_FAS * fas;
26421d9b32SPeter Brune   PetscErrorCode ierr;
27421d9b32SPeter Brune 
28421d9b32SPeter Brune   PetscFunctionBegin;
29421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
30421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
31421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
32421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
33421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
34421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
35421d9b32SPeter Brune 
36ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
37ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
38ed020824SBarry Smith 
39421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
40421d9b32SPeter Brune   snes->data                = (void*) fas;
41421d9b32SPeter Brune   fas->level                = 0;
42293a7e31SPeter Brune   fas->levels               = 1;
43ee78dd50SPeter Brune   fas->n_cycles             = 1;
44ee78dd50SPeter Brune   fas->max_up_it            = 1;
45ee78dd50SPeter Brune   fas->max_down_it          = 1;
46ee78dd50SPeter Brune   fas->upsmooth             = PETSC_NULL;
47ee78dd50SPeter Brune   fas->downsmooth           = PETSC_NULL;
48421d9b32SPeter Brune   fas->next                 = PETSC_NULL;
49421d9b32SPeter Brune   fas->interpolate          = PETSC_NULL;
50421d9b32SPeter Brune   fas->restrct              = PETSC_NULL;
51421d9b32SPeter Brune 
52421d9b32SPeter Brune   PetscFunctionReturn(0);
53421d9b32SPeter Brune }
54421d9b32SPeter Brune EXTERN_C_END
55421d9b32SPeter Brune 
56421d9b32SPeter Brune #undef __FUNCT__
57421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
58421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
59421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
60421d9b32SPeter Brune   PetscFunctionBegin;
61ee1fd11aSPeter Brune   *levels = fas->levels;
62421d9b32SPeter Brune   PetscFunctionReturn(0);
63421d9b32SPeter Brune }
64421d9b32SPeter Brune 
65421d9b32SPeter Brune #undef __FUNCT__
66421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
67421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
68421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
69421d9b32SPeter Brune   PetscInt levels = fas->level;
70421d9b32SPeter Brune   PetscInt i;
71421d9b32SPeter Brune   PetscFunctionBegin;
72421d9b32SPeter Brune   *lsnes = snes;
73421d9b32SPeter Brune   if (fas->level < level) {
74421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
75421d9b32SPeter Brune   }
76421d9b32SPeter Brune   if (level > levels - 1) {
77421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
78421d9b32SPeter Brune   }
79421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
80421d9b32SPeter Brune     *lsnes = fas->next;
81421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
82421d9b32SPeter Brune   }
83421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
84421d9b32SPeter Brune   PetscFunctionReturn(0);
85421d9b32SPeter Brune }
86421d9b32SPeter Brune 
87421d9b32SPeter Brune #undef __FUNCT__
88421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
89421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
90421d9b32SPeter Brune   PetscErrorCode ierr;
91421d9b32SPeter Brune   PetscInt i;
92421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
93421d9b32SPeter Brune   MPI_Comm comm;
94421d9b32SPeter Brune   PetscFunctionBegin;
95ee1fd11aSPeter Brune   comm = ((PetscObject)snes)->comm;
96ee1fd11aSPeter Brune   if (levels == fas->levels) {
97ee1fd11aSPeter Brune     if (!comms)
98ee1fd11aSPeter Brune       PetscFunctionReturn(0);
99ee1fd11aSPeter Brune   }
100421d9b32SPeter Brune   /* user has changed the number of levels; reset */
101421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
102421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
103421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
104ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
105421d9b32SPeter Brune   /* setup the finest level */
106421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
107421d9b32SPeter Brune     if (comms) comm = comms[i];
108421d9b32SPeter Brune     fas->level = i;
109421d9b32SPeter Brune     fas->levels = levels;
110ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
111421d9b32SPeter Brune     if (i > 0) {
112421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
113293a7e31SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
114421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
115421d9b32SPeter Brune       fas = (SNES_FAS *)fas->next->data;
116421d9b32SPeter Brune     }
117421d9b32SPeter Brune   }
118421d9b32SPeter Brune   PetscFunctionReturn(0);
119421d9b32SPeter Brune }
120421d9b32SPeter Brune 
121421d9b32SPeter Brune #undef __FUNCT__
122421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
123421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
124421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
125d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
126421d9b32SPeter Brune 
127421d9b32SPeter Brune   PetscFunctionBegin;
128421d9b32SPeter Brune   if (level > top_level)
129421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
130421d9b32SPeter Brune   /* get to the correct level */
131d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
132421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
133421d9b32SPeter Brune   }
134421d9b32SPeter Brune   if (fas->level != level)
135421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
136421d9b32SPeter Brune   fas->interpolate = mat;
137421d9b32SPeter Brune   PetscFunctionReturn(0);
138421d9b32SPeter Brune }
139421d9b32SPeter Brune 
140421d9b32SPeter Brune #undef __FUNCT__
141421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
142421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
143421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
144d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
145421d9b32SPeter Brune 
146421d9b32SPeter Brune   PetscFunctionBegin;
147421d9b32SPeter Brune   if (level > top_level)
148421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
149421d9b32SPeter Brune   /* get to the correct level */
150d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
151421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
152421d9b32SPeter Brune   }
153421d9b32SPeter Brune   if (fas->level != level)
154421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
155421d9b32SPeter Brune   fas->restrct = mat;
156421d9b32SPeter Brune   PetscFunctionReturn(0);
157421d9b32SPeter Brune }
158421d9b32SPeter Brune 
159421d9b32SPeter Brune #undef __FUNCT__
160421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
161421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
162421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
163d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
164421d9b32SPeter Brune 
165421d9b32SPeter Brune   PetscFunctionBegin;
166421d9b32SPeter Brune   if (level > top_level)
167421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
168421d9b32SPeter Brune   /* get to the correct level */
169d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
170421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
171421d9b32SPeter Brune   }
172421d9b32SPeter Brune   if (fas->level != level)
173421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
174421d9b32SPeter Brune   fas->rscale = rscale;
175421d9b32SPeter Brune   PetscFunctionReturn(0);
176421d9b32SPeter Brune }
177421d9b32SPeter Brune 
178421d9b32SPeter Brune #undef __FUNCT__
179421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
180421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
181421d9b32SPeter Brune {
18277df8cc4SPeter Brune   PetscErrorCode ierr = 0;
183421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
184421d9b32SPeter Brune 
185421d9b32SPeter Brune   PetscFunctionBegin;
186421d9b32SPeter Brune   /* destroy local data created in SNESSetup_FAS */
187ee78dd50SPeter Brune   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
188ee78dd50SPeter Brune   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
189293a7e31SPeter Brune   if (fas->interpolate == fas->restrct) {
190293a7e31SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
191293a7e31SPeter Brune     fas->restrct = PETSC_NULL;
192293a7e31SPeter Brune   } else {
193421d9b32SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
194421d9b32SPeter Brune     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
195293a7e31SPeter Brune   }
196421d9b32SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
197293a7e31SPeter Brune   /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */
198ee1fd11aSPeter Brune   if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr);
199fa9694d7SPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
200421d9b32SPeter Brune   PetscFunctionReturn(0);
201421d9b32SPeter Brune }
202421d9b32SPeter Brune 
203421d9b32SPeter Brune #undef __FUNCT__
204421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
205421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
206421d9b32SPeter Brune {
207421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
208421d9b32SPeter Brune   PetscErrorCode ierr;
209421d9b32SPeter Brune 
210421d9b32SPeter Brune   PetscFunctionBegin;
211421d9b32SPeter Brune   /* recursively resets and then destroys */
212421d9b32SPeter Brune   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
213421d9b32SPeter Brune   if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
214421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
215ee1fd11aSPeter Brune 
216421d9b32SPeter Brune   PetscFunctionReturn(0);
217421d9b32SPeter Brune }
218421d9b32SPeter Brune 
219421d9b32SPeter Brune #undef __FUNCT__
220421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
221421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
222421d9b32SPeter Brune {
223421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
224421d9b32SPeter Brune   PetscErrorCode ierr;
225421d9b32SPeter Brune 
226421d9b32SPeter Brune   PetscFunctionBegin;
227293a7e31SPeter Brune   ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
228421d9b32SPeter Brune   /* gets the solver ready for solution */
229421d9b32SPeter Brune   if (snes->dm) {
230421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
231421d9b32SPeter Brune     if (fas->next) {
232421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
233ee78dd50SPeter Brune         if (!fas->next->dm) {
234ee78dd50SPeter Brune           ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
235ee78dd50SPeter Brune           ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
236ee78dd50SPeter Brune         }
237421d9b32SPeter Brune         /* set the interpolation and restriction from the DM */
238ee78dd50SPeter Brune         if (!fas->interpolate) {
239421d9b32SPeter Brune           ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
240421d9b32SPeter Brune           fas->restrct = fas->interpolate;
241421d9b32SPeter Brune         }
242421d9b32SPeter Brune     }
243ee78dd50SPeter Brune     /* set the DMs of the pre and post-smoothers here */
244ee78dd50SPeter Brune     if (fas->upsmooth)  SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);
245ee78dd50SPeter Brune     if (fas->downsmooth)SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);
246421d9b32SPeter Brune   }
247fa9694d7SPeter Brune   if (fas->next) {
248fa9694d7SPeter Brune     /* gotta set up the solution vector for this to work */
249fa9694d7SPeter Brune     ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);
250fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
251fa9694d7SPeter Brune   }
252fa9694d7SPeter Brune   /* got to set them all up at once */
253421d9b32SPeter Brune   PetscFunctionReturn(0);
254421d9b32SPeter Brune }
255421d9b32SPeter Brune 
256421d9b32SPeter Brune #undef __FUNCT__
257421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
258421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
259421d9b32SPeter Brune {
260ee78dd50SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
261ee78dd50SPeter Brune   PetscInt levels = 1;
262ee78dd50SPeter Brune   PetscBool flg, monflg;
263421d9b32SPeter Brune   PetscErrorCode ierr;
264ee78dd50SPeter Brune   const char * def_smooth = SNESNRICHARDSON;
265ee78dd50SPeter Brune   char pre_type[256];
266ee78dd50SPeter Brune   char post_type[256];
267ee78dd50SPeter Brune   char                    monfilename[PETSC_MAX_PATH_LEN];
268421d9b32SPeter Brune 
269421d9b32SPeter Brune   PetscFunctionBegin;
270c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
271ee78dd50SPeter Brune 
272ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
273ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
274ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
275*c732cbdbSBarry Smith     if (!flg && snes->dm) {
276*c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
277*c732cbdbSBarry Smith       levels++;
278*c732cbdbSBarry Smith     }
279ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
280ee78dd50SPeter Brune   }
281ee78dd50SPeter Brune 
282ee78dd50SPeter Brune   /* type of pre and/or post smoothers -- set both at once */
283ee78dd50SPeter Brune   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
284ee78dd50SPeter Brune   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
285ee78dd50SPeter Brune   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
286ee78dd50SPeter Brune   if (flg) {
287ee78dd50SPeter Brune     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
288ee78dd50SPeter Brune   } else {
289ee78dd50SPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
290ee78dd50SPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
291ee78dd50SPeter Brune   }
292ee78dd50SPeter Brune 
293ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
294ee78dd50SPeter Brune   ierr = PetscOptionsInt("-snes_fas_up_it","Number of upsmoother iterations","PCMGSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
295ee78dd50SPeter Brune   ierr = PetscOptionsInt("-snes_fas_down_it","Number of downsmoother iterations","PCMGSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
296ee78dd50SPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","PCMGSetNumberSmoothUp",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
297ee78dd50SPeter Brune 
298ee78dd50SPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor for smoothers","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
299ee78dd50SPeter Brune 
300ee78dd50SPeter Brune   /* other options for the coarsest level */
301ee78dd50SPeter Brune   if (fas->level == 0) {
302ee78dd50SPeter Brune     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
303ee78dd50SPeter Brune     if (flg) {
304ee78dd50SPeter Brune       ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
305ee78dd50SPeter Brune     } else {
306ee78dd50SPeter Brune       ierr = PetscOptionsList("-snes_fas_coarse_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
307ee78dd50SPeter Brune       ierr = PetscOptionsList("-snes_fas_coarse_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
308ee78dd50SPeter Brune     }
309ee78dd50SPeter Brune   }
310ee78dd50SPeter Brune 
311421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
312ee78dd50SPeter Brune   /* setup from the determined types if the smoothers don't exist */
313ee78dd50SPeter Brune   if (!fas->upsmooth) {
31467339d5cSBarry Smith     const char     *prefix;
31567339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
316ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
31767339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
31867339d5cSBarry Smith     if (fas->level != 0) {
31967339d5cSBarry Smith       ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_");CHKERRQ(ierr);
32067339d5cSBarry Smith     } else {
32167339d5cSBarry Smith       ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_coarse_");CHKERRQ(ierr);
32267339d5cSBarry Smith     }
323293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
324ee78dd50SPeter Brune     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
325ee78dd50SPeter Brune   }
326293a7e31SPeter Brune   if (fas->upsmooth) ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
327293a7e31SPeter Brune   if (!fas->downsmooth && fas->level != 0) {
32867339d5cSBarry Smith     const char     *prefix;
32967339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
330ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
33167339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
33267339d5cSBarry Smith     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_");CHKERRQ(ierr);
333293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
334ee78dd50SPeter Brune     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
335ee78dd50SPeter Brune   }
336293a7e31SPeter Brune   if (fas->downsmooth) ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
337ee78dd50SPeter Brune 
338ee78dd50SPeter Brune   if (monflg) {
339293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
340293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
341ee78dd50SPeter Brune   }
342ee78dd50SPeter Brune 
343ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
344ee1fd11aSPeter Brune   if (fas->next)ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);
345421d9b32SPeter Brune   PetscFunctionReturn(0);
346421d9b32SPeter Brune }
347421d9b32SPeter Brune 
348421d9b32SPeter Brune #undef __FUNCT__
349421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
350421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
351421d9b32SPeter Brune {
352421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
353421d9b32SPeter Brune   PetscBool      iascii;
354421d9b32SPeter Brune   PetscErrorCode ierr;
355421d9b32SPeter Brune   PetscInt levels = fas->levels;
356421d9b32SPeter Brune   PetscInt i;
357421d9b32SPeter Brune 
358421d9b32SPeter Brune   PetscFunctionBegin;
359421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
360421d9b32SPeter Brune   if (iascii) {
361421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
362421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
363421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
364421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
365ee78dd50SPeter Brune       if (fas->upsmooth) {
366ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
367421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
368ee78dd50SPeter Brune         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
369421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
370421d9b32SPeter Brune       } else {
371ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
372421d9b32SPeter Brune       }
373ee78dd50SPeter Brune       if (fas->downsmooth) {
374ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
375421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
376ee78dd50SPeter Brune         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
377421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
378421d9b32SPeter Brune       } else {
379ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
380421d9b32SPeter Brune       }
381421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
382421d9b32SPeter Brune     }
383421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
384421d9b32SPeter Brune   } else {
385421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
386421d9b32SPeter Brune   }
387421d9b32SPeter Brune   PetscFunctionReturn(0);
388421d9b32SPeter Brune }
389421d9b32SPeter Brune 
390421d9b32SPeter Brune /*
391fa9694d7SPeter Brune 
392fa9694d7SPeter Brune Defines the FAS cycle as:
393fa9694d7SPeter Brune 
394fa9694d7SPeter Brune fine problem: F(x) = 0
395fa9694d7SPeter Brune coarse problem: F^c(x) = b^c
396fa9694d7SPeter Brune 
397fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
398fa9694d7SPeter Brune 
399fa9694d7SPeter Brune correction:
400fa9694d7SPeter Brune 
401fa9694d7SPeter Brune x = x + I(x^c - Rx)
402fa9694d7SPeter Brune 
403421d9b32SPeter Brune  */
404fa9694d7SPeter Brune 
405421d9b32SPeter Brune #undef __FUNCT__
406421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private"
407c90fad12SPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec B, Vec X, Vec F) {
408fa9694d7SPeter Brune 
409fa9694d7SPeter Brune   PetscErrorCode ierr;
410fa9694d7SPeter Brune   Vec X_c, Xo_c, F_c, B_c;
411fa9694d7SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
412293a7e31SPeter Brune   PetscInt i;
413421d9b32SPeter Brune 
414421d9b32SPeter Brune   PetscFunctionBegin;
415fa9694d7SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
416ee78dd50SPeter Brune   if (fas->upsmooth) {
417ee78dd50SPeter Brune     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
418c90fad12SPeter Brune   } else if (snes->pc) {
419c90fad12SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
420fe6f9142SPeter Brune   }
421fa9694d7SPeter Brune   if (fas->next) {
422293a7e31SPeter Brune     for (i = 0; i < fas->n_cycles; i++) {
423ee78dd50SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
424c90fad12SPeter Brune       X_c  = fas->next->vec_sol;
425293a7e31SPeter Brune       Xo_c = fas->next->work[0];
426c90fad12SPeter Brune       F_c  = fas->next->vec_func;
427293a7e31SPeter Brune       B_c  = fas->next->work[1];
428293a7e31SPeter Brune       /* inject the solution to coarse */
429c90fad12SPeter Brune       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
430fa9694d7SPeter Brune       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
431c90fad12SPeter Brune       if (B) {
432293a7e31SPeter Brune         ierr = VecAYPX(F, -1.0, B);CHKERRQ(ierr); /* F = B - F (defect) */
433c90fad12SPeter Brune       } else {
434293a7e31SPeter Brune         ierr = VecScale(F, -1.0);CHKERRQ(ierr);
435c90fad12SPeter Brune       }
436293a7e31SPeter Brune 
437293a7e31SPeter Brune       /* restrict the defect */
438293a7e31SPeter Brune       ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
439293a7e31SPeter Brune       ierr = VecPointwiseMult(B_c,  fas->rscale, B_c);CHKERRQ(ierr);
440293a7e31SPeter Brune 
441c90fad12SPeter Brune       /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
442ee78dd50SPeter Brune       fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
443c90fad12SPeter Brune       ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
444293a7e31SPeter Brune       ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
445c90fad12SPeter Brune 
446ee78dd50SPeter Brune       /* set initial guess of the coarse problem to the projected fine solution */
447ee78dd50SPeter Brune       ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
448c90fad12SPeter Brune 
449c90fad12SPeter Brune       /* recurse to the next level */
450c90fad12SPeter Brune       ierr = FASCycle_Private(fas->next, B_c, X_c, F_c);CHKERRQ(ierr);
451ee78dd50SPeter Brune 
452fa9694d7SPeter Brune       /* correct as x <- x + I(x^c - Rx)*/
453fa9694d7SPeter Brune       ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
454ee78dd50SPeter Brune       ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr);
455293a7e31SPeter Brune     }
456fa9694d7SPeter Brune   }
457ee78dd50SPeter Brune     /* down-smooth -- just update using the down-smoother */
458c90fad12SPeter Brune   if (fas->level != 0) {
459ee78dd50SPeter Brune     if (fas->downsmooth) {
460ee78dd50SPeter Brune       ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
461c90fad12SPeter Brune     } else if (snes->pc) {
462c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
463c90fad12SPeter Brune     }
464fe6f9142SPeter Brune   }
465fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
466fa9694d7SPeter Brune 
467fa9694d7SPeter Brune   PetscFunctionReturn(0);
468421d9b32SPeter Brune }
469421d9b32SPeter Brune 
470421d9b32SPeter Brune #undef __FUNCT__
471c90fad12SPeter Brune #define __FUNCT__ "FASInitialGuess_Private"
472293a7e31SPeter Brune PetscErrorCode FASInitialGuess_Private(SNES snes, Vec B, Vec X) {
473c90fad12SPeter Brune 
474293a7e31SPeter Brune   PetscErrorCode ierr;
475293a7e31SPeter Brune   Vec X_c, B_c;
476ee1fd11aSPeter Brune   SNES_FAS * fas;
477293a7e31SPeter Brune 
478293a7e31SPeter Brune   PetscFunctionBegin;
479ee1fd11aSPeter Brune   fas = (SNES_FAS *)snes->data;
480293a7e31SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
481293a7e31SPeter Brune   if (fas->level == 0) {
482293a7e31SPeter Brune     if (fas->upsmooth) {
483293a7e31SPeter Brune       ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
484293a7e31SPeter Brune     } else if (snes->pc) {
485293a7e31SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
486293a7e31SPeter Brune     }
487293a7e31SPeter Brune   }
488293a7e31SPeter Brune   if (fas->next) {
489293a7e31SPeter Brune     X_c  = fas->next->vec_sol;
490293a7e31SPeter Brune     B_c  = fas->next->work[0];
491293a7e31SPeter Brune     /* inject the solution to coarse */
492293a7e31SPeter Brune     ierr = MatRestrict(fas->restrct, X, X_c);CHKERRQ(ierr);
493293a7e31SPeter Brune     ierr = VecPointwiseMult(X_c, fas->rscale, X_c);CHKERRQ(ierr);
494293a7e31SPeter Brune     if (B) {
495293a7e31SPeter Brune       ierr = MatRestrict(fas->restrct, B, B_c);CHKERRQ(ierr);
496293a7e31SPeter Brune       ierr = VecPointwiseMult(B_c, fas->rscale, B_c);CHKERRQ(ierr);
497293a7e31SPeter Brune     } else {
498293a7e31SPeter Brune       B_c = PETSC_NULL;
499293a7e31SPeter Brune     }
500293a7e31SPeter Brune     /* recurse to the next level */
501293a7e31SPeter Brune     ierr = FASInitialGuess_Private(fas->next, B_c, X_c);CHKERRQ(ierr);
502293a7e31SPeter Brune     ierr = MatInterpolate(fas->interpolate, X_c, X);CHKERRQ(ierr);
503293a7e31SPeter Brune   }
504293a7e31SPeter Brune   /* down-smooth -- just update using the down-smoother */
505293a7e31SPeter Brune   if (fas->level != 0) {
506293a7e31SPeter Brune     if (fas->downsmooth) {
507293a7e31SPeter Brune       ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
508293a7e31SPeter Brune     } else if (snes->pc) {
509293a7e31SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
510293a7e31SPeter Brune     }
511293a7e31SPeter Brune   }
512293a7e31SPeter Brune   PetscFunctionReturn(0);
513293a7e31SPeter Brune }
514c90fad12SPeter Brune 
515c90fad12SPeter Brune #undef __FUNCT__
516421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
517421d9b32SPeter Brune 
518421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
519421d9b32SPeter Brune {
520fa9694d7SPeter Brune   PetscErrorCode ierr;
521fe6f9142SPeter Brune   PetscInt i, maxits;
522fe6f9142SPeter Brune   Vec X, F, B;
523fe6f9142SPeter Brune   PetscReal fnorm;
524421d9b32SPeter Brune   PetscFunctionBegin;
525fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
526fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
527fa9694d7SPeter Brune   X = snes->vec_sol;
528fa9694d7SPeter Brune   F = snes->vec_func;
529fe6f9142SPeter Brune   B = snes->vec_rhs;
530293a7e31SPeter Brune 
531293a7e31SPeter Brune   /*norm setup */
532fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
533fe6f9142SPeter Brune   snes->iter = 0;
534fe6f9142SPeter Brune   snes->norm = 0.;
535fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
536fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
537fe6f9142SPeter Brune   if (snes->domainerror) {
538fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
539fe6f9142SPeter Brune     PetscFunctionReturn(0);
540fe6f9142SPeter Brune   }
541fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
542fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
543fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
544fe6f9142SPeter Brune   snes->norm = fnorm;
545fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
546fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
547fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
548fe6f9142SPeter Brune 
549fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
550fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
551fe6f9142SPeter Brune   /* test convergence */
552fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
553fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
554fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
555fe6f9142SPeter Brune     /* Call general purpose update function */
556fe6f9142SPeter Brune     if (snes->ops->update) {
557fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
558fe6f9142SPeter Brune     }
559c90fad12SPeter Brune     ierr = FASCycle_Private(snes, B, X, F);CHKERRQ(ierr);
560c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
561c90fad12SPeter Brune     /* Monitor convergence */
562c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
563c90fad12SPeter Brune     snes->iter = i+1;
564c90fad12SPeter Brune     snes->norm = fnorm;
565c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
566c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
567c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
568c90fad12SPeter Brune     /* Test for convergence */
569c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
570c90fad12SPeter Brune     if (snes->reason) break;
571fe6f9142SPeter Brune   }
572fe6f9142SPeter Brune   if (i == maxits) {
573fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
574fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
575fe6f9142SPeter Brune   }
576421d9b32SPeter Brune   PetscFunctionReturn(0);
577421d9b32SPeter Brune }
578