xref: /petsc/src/snes/impls/fas/fas.c (revision ee1fd11aa01d29bf27f6b93858a93914927104a5)
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 
36421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
37421d9b32SPeter Brune   snes->data                = (void*) fas;
38421d9b32SPeter Brune   fas->level                = 0;
39293a7e31SPeter Brune   fas->levels               = 1;
40ee78dd50SPeter Brune   fas->n_cycles             = 1;
41ee78dd50SPeter Brune   fas->max_up_it            = 1;
42ee78dd50SPeter Brune   fas->max_down_it          = 1;
43ee78dd50SPeter Brune   fas->upsmooth             = PETSC_NULL;
44ee78dd50SPeter Brune   fas->downsmooth           = PETSC_NULL;
45421d9b32SPeter Brune   fas->next                 = PETSC_NULL;
46421d9b32SPeter Brune   fas->interpolate          = PETSC_NULL;
47421d9b32SPeter Brune   fas->restrct              = PETSC_NULL;
48421d9b32SPeter Brune 
49421d9b32SPeter Brune   PetscFunctionReturn(0);
50421d9b32SPeter Brune }
51421d9b32SPeter Brune EXTERN_C_END
52421d9b32SPeter Brune 
53421d9b32SPeter Brune #undef __FUNCT__
54421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
55421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
56421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
57421d9b32SPeter Brune   PetscFunctionBegin;
58*ee1fd11aSPeter Brune   *levels = fas->levels;
59421d9b32SPeter Brune   PetscFunctionReturn(0);
60421d9b32SPeter Brune }
61421d9b32SPeter Brune 
62421d9b32SPeter Brune #undef __FUNCT__
63421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
64421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
65421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
66421d9b32SPeter Brune   PetscInt levels = fas->level;
67421d9b32SPeter Brune   PetscInt i;
68421d9b32SPeter Brune   PetscFunctionBegin;
69421d9b32SPeter Brune   *lsnes = snes;
70421d9b32SPeter Brune   if (fas->level < level) {
71421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
72421d9b32SPeter Brune   }
73421d9b32SPeter Brune   if (level > levels - 1) {
74421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
75421d9b32SPeter Brune   }
76421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
77421d9b32SPeter Brune     *lsnes = fas->next;
78421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
79421d9b32SPeter Brune   }
80421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
81421d9b32SPeter Brune   PetscFunctionReturn(0);
82421d9b32SPeter Brune }
83421d9b32SPeter Brune 
84421d9b32SPeter Brune #undef __FUNCT__
85421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
86421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
87421d9b32SPeter Brune   PetscErrorCode ierr;
88421d9b32SPeter Brune   PetscInt i;
89421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
90421d9b32SPeter Brune   MPI_Comm comm;
91421d9b32SPeter Brune   PetscFunctionBegin;
92*ee1fd11aSPeter Brune   comm = ((PetscObject)snes)->comm;
93*ee1fd11aSPeter Brune   if (levels == fas->levels) {
94*ee1fd11aSPeter Brune     if (!comms)
95*ee1fd11aSPeter Brune       PetscFunctionReturn(0);
96*ee1fd11aSPeter Brune   }
97421d9b32SPeter Brune   /* user has changed the number of levels; reset */
98421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
99421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
100421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
101*ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
102421d9b32SPeter Brune   /* setup the finest level */
103421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
104421d9b32SPeter Brune     if (comms) comm = comms[i];
105421d9b32SPeter Brune     fas->level = i;
106421d9b32SPeter Brune     fas->levels = levels;
107*ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
108421d9b32SPeter Brune     if (i > 0) {
109421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
110293a7e31SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
111421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
112421d9b32SPeter Brune       fas = (SNES_FAS *)fas->next->data;
113421d9b32SPeter Brune     }
114421d9b32SPeter Brune   }
115421d9b32SPeter Brune   PetscFunctionReturn(0);
116421d9b32SPeter Brune }
117421d9b32SPeter Brune 
118421d9b32SPeter Brune #undef __FUNCT__
119421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
120421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
121421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
122d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
123421d9b32SPeter Brune 
124421d9b32SPeter Brune   PetscFunctionBegin;
125421d9b32SPeter Brune   if (level > top_level)
126421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
127421d9b32SPeter Brune   /* get to the correct level */
128d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
129421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
130421d9b32SPeter Brune   }
131421d9b32SPeter Brune   if (fas->level != level)
132421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
133421d9b32SPeter Brune   fas->interpolate = mat;
134421d9b32SPeter Brune   PetscFunctionReturn(0);
135421d9b32SPeter Brune }
136421d9b32SPeter Brune 
137421d9b32SPeter Brune #undef __FUNCT__
138421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
139421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
140421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
141d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
142421d9b32SPeter Brune 
143421d9b32SPeter Brune   PetscFunctionBegin;
144421d9b32SPeter Brune   if (level > top_level)
145421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
146421d9b32SPeter Brune   /* get to the correct level */
147d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
148421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
149421d9b32SPeter Brune   }
150421d9b32SPeter Brune   if (fas->level != level)
151421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
152421d9b32SPeter Brune   fas->restrct = mat;
153421d9b32SPeter Brune   PetscFunctionReturn(0);
154421d9b32SPeter Brune }
155421d9b32SPeter Brune 
156421d9b32SPeter Brune #undef __FUNCT__
157421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
158421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
159421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
160d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
161421d9b32SPeter Brune 
162421d9b32SPeter Brune   PetscFunctionBegin;
163421d9b32SPeter Brune   if (level > top_level)
164421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
165421d9b32SPeter Brune   /* get to the correct level */
166d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
167421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
168421d9b32SPeter Brune   }
169421d9b32SPeter Brune   if (fas->level != level)
170421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
171421d9b32SPeter Brune   fas->rscale = rscale;
172421d9b32SPeter Brune   PetscFunctionReturn(0);
173421d9b32SPeter Brune }
174421d9b32SPeter Brune 
175421d9b32SPeter Brune #undef __FUNCT__
176421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
177421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
178421d9b32SPeter Brune {
17977df8cc4SPeter Brune   PetscErrorCode ierr = 0;
180421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
181421d9b32SPeter Brune 
182421d9b32SPeter Brune   PetscFunctionBegin;
183421d9b32SPeter Brune   /* destroy local data created in SNESSetup_FAS */
184ee78dd50SPeter Brune   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
185ee78dd50SPeter Brune   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
186293a7e31SPeter Brune   if (fas->interpolate == fas->restrct) {
187293a7e31SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
188293a7e31SPeter Brune     fas->restrct = PETSC_NULL;
189293a7e31SPeter Brune   } else {
190421d9b32SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
191421d9b32SPeter Brune     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
192293a7e31SPeter Brune   }
193421d9b32SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
194293a7e31SPeter Brune   /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */
195*ee1fd11aSPeter Brune   if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr);
196fa9694d7SPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
197421d9b32SPeter Brune   PetscFunctionReturn(0);
198421d9b32SPeter Brune }
199421d9b32SPeter Brune 
200421d9b32SPeter Brune #undef __FUNCT__
201421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
202421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
203421d9b32SPeter Brune {
204421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
205421d9b32SPeter Brune   PetscErrorCode ierr;
206421d9b32SPeter Brune 
207421d9b32SPeter Brune   PetscFunctionBegin;
208421d9b32SPeter Brune   /* recursively resets and then destroys */
209421d9b32SPeter Brune   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
210421d9b32SPeter Brune   if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
211421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
212*ee1fd11aSPeter Brune 
213421d9b32SPeter Brune   PetscFunctionReturn(0);
214421d9b32SPeter Brune }
215421d9b32SPeter Brune 
216421d9b32SPeter Brune #undef __FUNCT__
217421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
218421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
219421d9b32SPeter Brune {
220421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
221421d9b32SPeter Brune   PetscErrorCode ierr;
222421d9b32SPeter Brune 
223421d9b32SPeter Brune   PetscFunctionBegin;
224293a7e31SPeter Brune   ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
225421d9b32SPeter Brune   /* gets the solver ready for solution */
226421d9b32SPeter Brune   if (snes->dm) {
227421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
228421d9b32SPeter Brune     if (fas->next) {
229421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
230ee78dd50SPeter Brune         if (!fas->next->dm) {
231ee78dd50SPeter Brune           ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
232ee78dd50SPeter Brune           ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
233ee78dd50SPeter Brune         }
234421d9b32SPeter Brune         /* set the interpolation and restriction from the DM */
235ee78dd50SPeter Brune         if (!fas->interpolate) {
236421d9b32SPeter Brune           ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
237421d9b32SPeter Brune           fas->restrct = fas->interpolate;
238421d9b32SPeter Brune         }
239421d9b32SPeter Brune     }
240ee78dd50SPeter Brune     /* set the DMs of the pre and post-smoothers here */
241ee78dd50SPeter Brune     if (fas->upsmooth)  SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);
242ee78dd50SPeter Brune     if (fas->downsmooth)SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);
243421d9b32SPeter Brune   }
244fa9694d7SPeter Brune   if (fas->next) {
245fa9694d7SPeter Brune     /* gotta set up the solution vector for this to work */
246fa9694d7SPeter Brune     ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);
247fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
248fa9694d7SPeter Brune   }
249fa9694d7SPeter Brune   /* got to set them all up at once */
250421d9b32SPeter Brune   PetscFunctionReturn(0);
251421d9b32SPeter Brune }
252421d9b32SPeter Brune 
253421d9b32SPeter Brune #undef __FUNCT__
254421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
255421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
256421d9b32SPeter Brune {
257ee78dd50SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
258ee78dd50SPeter Brune   PetscInt levels = 1;
259ee78dd50SPeter Brune   PetscBool flg, monflg;
260421d9b32SPeter Brune   PetscErrorCode ierr;
261ee78dd50SPeter Brune   const char * def_smooth = SNESNRICHARDSON;
262ee78dd50SPeter Brune   char pre_type[256];
263ee78dd50SPeter Brune   char post_type[256];
264ee78dd50SPeter Brune   char                    monfilename[PETSC_MAX_PATH_LEN];
265421d9b32SPeter Brune 
266421d9b32SPeter Brune   PetscFunctionBegin;
267c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
268ee78dd50SPeter Brune 
269ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
270ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
271ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
272ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
273ee78dd50SPeter Brune   }
274ee78dd50SPeter Brune 
275ee78dd50SPeter Brune   /* type of pre and/or post smoothers -- set both at once */
276ee78dd50SPeter Brune   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
277ee78dd50SPeter Brune   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
278ee78dd50SPeter Brune   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
279ee78dd50SPeter Brune   if (flg) {
280ee78dd50SPeter Brune     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
281ee78dd50SPeter Brune   } else {
282ee78dd50SPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
283ee78dd50SPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
284ee78dd50SPeter Brune   }
285ee78dd50SPeter Brune 
286ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
287ee78dd50SPeter Brune   ierr = PetscOptionsInt("-snes_fas_up_it","Number of upsmoother iterations","PCMGSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
288ee78dd50SPeter Brune   ierr = PetscOptionsInt("-snes_fas_down_it","Number of downsmoother iterations","PCMGSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
289ee78dd50SPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","PCMGSetNumberSmoothUp",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
290ee78dd50SPeter Brune 
291ee78dd50SPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor for smoothers","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
292ee78dd50SPeter Brune 
293ee78dd50SPeter Brune   /* other options for the coarsest level */
294ee78dd50SPeter Brune   if (fas->level == 0) {
295ee78dd50SPeter Brune     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
296ee78dd50SPeter Brune     if (flg) {
297ee78dd50SPeter Brune       ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
298ee78dd50SPeter Brune     } else {
299ee78dd50SPeter Brune       ierr = PetscOptionsList("-snes_fas_coarse_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
300ee78dd50SPeter Brune       ierr = PetscOptionsList("-snes_fas_coarse_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
301ee78dd50SPeter Brune     }
302ee78dd50SPeter Brune   }
303ee78dd50SPeter Brune 
304421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
305ee78dd50SPeter Brune   /* setup from the determined types if the smoothers don't exist */
306ee78dd50SPeter Brune   if (!fas->upsmooth) {
307ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
308293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
309ee78dd50SPeter Brune     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
310ee78dd50SPeter Brune   }
311293a7e31SPeter Brune   if (fas->upsmooth) ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
312293a7e31SPeter Brune   if (!fas->downsmooth && fas->level != 0) {
313ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
314293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
315ee78dd50SPeter Brune    ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
316ee78dd50SPeter Brune   }
317293a7e31SPeter Brune   if (fas->downsmooth) ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
318ee78dd50SPeter Brune 
319ee78dd50SPeter Brune   if (monflg) {
320293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
321293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
322ee78dd50SPeter Brune   }
323ee78dd50SPeter Brune 
324ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
325*ee1fd11aSPeter Brune   if (fas->next)ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);
326421d9b32SPeter Brune   PetscFunctionReturn(0);
327421d9b32SPeter Brune }
328421d9b32SPeter Brune 
329421d9b32SPeter Brune #undef __FUNCT__
330421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
331421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
332421d9b32SPeter Brune {
333421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
334421d9b32SPeter Brune   PetscBool      iascii;
335421d9b32SPeter Brune   PetscErrorCode ierr;
336421d9b32SPeter Brune   PetscInt levels = fas->levels;
337421d9b32SPeter Brune   PetscInt i;
338421d9b32SPeter Brune 
339421d9b32SPeter Brune   PetscFunctionBegin;
340421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
341421d9b32SPeter Brune   if (iascii) {
342421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
343421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
344421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
345421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
346ee78dd50SPeter Brune       if (fas->upsmooth) {
347ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
348421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
349ee78dd50SPeter Brune         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
350421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
351421d9b32SPeter Brune       } else {
352ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
353421d9b32SPeter Brune       }
354ee78dd50SPeter Brune       if (fas->downsmooth) {
355ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
356421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
357ee78dd50SPeter Brune         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
358421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
359421d9b32SPeter Brune       } else {
360ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
361421d9b32SPeter Brune       }
362421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
363421d9b32SPeter Brune     }
364421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
365421d9b32SPeter Brune   } else {
366421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
367421d9b32SPeter Brune   }
368421d9b32SPeter Brune   PetscFunctionReturn(0);
369421d9b32SPeter Brune }
370421d9b32SPeter Brune 
371421d9b32SPeter Brune /*
372fa9694d7SPeter Brune 
373fa9694d7SPeter Brune Defines the FAS cycle as:
374fa9694d7SPeter Brune 
375fa9694d7SPeter Brune fine problem: F(x) = 0
376fa9694d7SPeter Brune coarse problem: F^c(x) = b^c
377fa9694d7SPeter Brune 
378fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
379fa9694d7SPeter Brune 
380fa9694d7SPeter Brune correction:
381fa9694d7SPeter Brune 
382fa9694d7SPeter Brune x = x + I(x^c - Rx)
383fa9694d7SPeter Brune 
384421d9b32SPeter Brune  */
385fa9694d7SPeter Brune 
386421d9b32SPeter Brune #undef __FUNCT__
387421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private"
388c90fad12SPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec B, Vec X, Vec F) {
389fa9694d7SPeter Brune 
390fa9694d7SPeter Brune   PetscErrorCode ierr;
391fa9694d7SPeter Brune   Vec X_c, Xo_c, F_c, B_c;
392fa9694d7SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
393293a7e31SPeter Brune   PetscInt i;
394421d9b32SPeter Brune 
395421d9b32SPeter Brune   PetscFunctionBegin;
396fa9694d7SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
397ee78dd50SPeter Brune   if (fas->upsmooth) {
398ee78dd50SPeter Brune     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
399c90fad12SPeter Brune   } else if (snes->pc) {
400c90fad12SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
401fe6f9142SPeter Brune   }
402fa9694d7SPeter Brune   if (fas->next) {
403293a7e31SPeter Brune     for (i = 0; i < fas->n_cycles; i++) {
404ee78dd50SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
405c90fad12SPeter Brune       X_c  = fas->next->vec_sol;
406293a7e31SPeter Brune       Xo_c = fas->next->work[0];
407c90fad12SPeter Brune       F_c  = fas->next->vec_func;
408293a7e31SPeter Brune       B_c  = fas->next->work[1];
409293a7e31SPeter Brune       /* inject the solution to coarse */
410c90fad12SPeter Brune       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
411fa9694d7SPeter Brune       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
412c90fad12SPeter Brune       if (B) {
413293a7e31SPeter Brune         ierr = VecAYPX(F, -1.0, B);CHKERRQ(ierr); /* F = B - F (defect) */
414c90fad12SPeter Brune       } else {
415293a7e31SPeter Brune         ierr = VecScale(F, -1.0);CHKERRQ(ierr);
416c90fad12SPeter Brune       }
417293a7e31SPeter Brune 
418293a7e31SPeter Brune       /* restrict the defect */
419293a7e31SPeter Brune       ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
420293a7e31SPeter Brune       ierr = VecPointwiseMult(B_c,  fas->rscale, B_c);CHKERRQ(ierr);
421293a7e31SPeter Brune 
422c90fad12SPeter Brune       /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
423ee78dd50SPeter Brune       fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
424c90fad12SPeter Brune       ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
425293a7e31SPeter Brune       ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
426c90fad12SPeter Brune 
427ee78dd50SPeter Brune       /* set initial guess of the coarse problem to the projected fine solution */
428ee78dd50SPeter Brune       ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
429c90fad12SPeter Brune 
430c90fad12SPeter Brune       /* recurse to the next level */
431c90fad12SPeter Brune       ierr = FASCycle_Private(fas->next, B_c, X_c, F_c);CHKERRQ(ierr);
432ee78dd50SPeter Brune 
433fa9694d7SPeter Brune       /* correct as x <- x + I(x^c - Rx)*/
434fa9694d7SPeter Brune       ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
435ee78dd50SPeter Brune       ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr);
436293a7e31SPeter Brune     }
437fa9694d7SPeter Brune   }
438ee78dd50SPeter Brune     /* down-smooth -- just update using the down-smoother */
439c90fad12SPeter Brune   if (fas->level != 0) {
440ee78dd50SPeter Brune     if (fas->downsmooth) {
441ee78dd50SPeter Brune       ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
442c90fad12SPeter Brune     } else if (snes->pc) {
443c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
444c90fad12SPeter Brune     }
445fe6f9142SPeter Brune   }
446fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
447fa9694d7SPeter Brune 
448fa9694d7SPeter Brune   PetscFunctionReturn(0);
449421d9b32SPeter Brune }
450421d9b32SPeter Brune 
451421d9b32SPeter Brune #undef __FUNCT__
452c90fad12SPeter Brune #define __FUNCT__ "FASInitialGuess_Private"
453293a7e31SPeter Brune PetscErrorCode FASInitialGuess_Private(SNES snes, Vec B, Vec X) {
454c90fad12SPeter Brune 
455293a7e31SPeter Brune   PetscErrorCode ierr;
456293a7e31SPeter Brune   Vec X_c, B_c;
457*ee1fd11aSPeter Brune   SNES_FAS * fas;
458293a7e31SPeter Brune 
459293a7e31SPeter Brune   PetscFunctionBegin;
460*ee1fd11aSPeter Brune   fas = (SNES_FAS *)snes->data;
461293a7e31SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
462293a7e31SPeter Brune   if (fas->level == 0) {
463293a7e31SPeter Brune     if (fas->upsmooth) {
464293a7e31SPeter Brune       ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
465293a7e31SPeter Brune     } else if (snes->pc) {
466293a7e31SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
467293a7e31SPeter Brune     }
468293a7e31SPeter Brune   }
469293a7e31SPeter Brune   if (fas->next) {
470293a7e31SPeter Brune     X_c  = fas->next->vec_sol;
471293a7e31SPeter Brune     B_c  = fas->next->work[0];
472293a7e31SPeter Brune     /* inject the solution to coarse */
473293a7e31SPeter Brune     ierr = MatRestrict(fas->restrct, X, X_c);CHKERRQ(ierr);
474293a7e31SPeter Brune     ierr = VecPointwiseMult(X_c, fas->rscale, X_c);CHKERRQ(ierr);
475293a7e31SPeter Brune     if (B) {
476293a7e31SPeter Brune       ierr = MatRestrict(fas->restrct, B, B_c);CHKERRQ(ierr);
477293a7e31SPeter Brune       ierr = VecPointwiseMult(B_c, fas->rscale, B_c);CHKERRQ(ierr);
478293a7e31SPeter Brune     } else {
479293a7e31SPeter Brune       B_c = PETSC_NULL;
480293a7e31SPeter Brune     }
481293a7e31SPeter Brune     /* recurse to the next level */
482293a7e31SPeter Brune     ierr = FASInitialGuess_Private(fas->next, B_c, X_c);CHKERRQ(ierr);
483293a7e31SPeter Brune     ierr = MatInterpolate(fas->interpolate, X_c, X);CHKERRQ(ierr);
484293a7e31SPeter Brune   }
485293a7e31SPeter Brune   /* down-smooth -- just update using the down-smoother */
486293a7e31SPeter Brune   if (fas->level != 0) {
487293a7e31SPeter Brune     if (fas->downsmooth) {
488293a7e31SPeter Brune       ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
489293a7e31SPeter Brune     } else if (snes->pc) {
490293a7e31SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
491293a7e31SPeter Brune     }
492293a7e31SPeter Brune   }
493293a7e31SPeter Brune   PetscFunctionReturn(0);
494293a7e31SPeter Brune }
495c90fad12SPeter Brune 
496c90fad12SPeter Brune #undef __FUNCT__
497421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
498421d9b32SPeter Brune 
499421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
500421d9b32SPeter Brune {
501fa9694d7SPeter Brune   PetscErrorCode ierr;
502fe6f9142SPeter Brune   PetscInt i, maxits;
503fe6f9142SPeter Brune   Vec X, F, B;
504fe6f9142SPeter Brune   PetscReal fnorm;
505421d9b32SPeter Brune   PetscFunctionBegin;
506fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
507fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
508fa9694d7SPeter Brune   X = snes->vec_sol;
509fa9694d7SPeter Brune   F = snes->vec_func;
510fe6f9142SPeter Brune   B = snes->vec_rhs;
511293a7e31SPeter Brune 
512293a7e31SPeter Brune   /*norm setup */
513fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
514fe6f9142SPeter Brune   snes->iter = 0;
515fe6f9142SPeter Brune   snes->norm = 0.;
516fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
517fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
518fe6f9142SPeter Brune   if (snes->domainerror) {
519fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
520fe6f9142SPeter Brune     PetscFunctionReturn(0);
521fe6f9142SPeter Brune   }
522fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
523fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
524fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
525fe6f9142SPeter Brune   snes->norm = fnorm;
526fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
527fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
528fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
529fe6f9142SPeter Brune 
530fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
531fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
532fe6f9142SPeter Brune   /* test convergence */
533fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
534fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
535fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
536fe6f9142SPeter Brune     /* Call general purpose update function */
537fe6f9142SPeter Brune     if (snes->ops->update) {
538fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
539fe6f9142SPeter Brune     }
540c90fad12SPeter Brune     ierr = FASCycle_Private(snes, B, X, F);CHKERRQ(ierr);
541c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
542c90fad12SPeter Brune     /* Monitor convergence */
543c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
544c90fad12SPeter Brune     snes->iter = i+1;
545c90fad12SPeter Brune     snes->norm = fnorm;
546c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
547c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
548c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
549c90fad12SPeter Brune     /* Test for convergence */
550c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
551c90fad12SPeter Brune     if (snes->reason) break;
552fe6f9142SPeter Brune   }
553fe6f9142SPeter Brune   if (i == maxits) {
554fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
555fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
556fe6f9142SPeter Brune   }
557421d9b32SPeter Brune   PetscFunctionReturn(0);
558421d9b32SPeter Brune }
559