xref: /petsc/src/snes/impls/fas/fas.c (revision 794bee337f0d5af71ccfe93f8f7f88922f7fb150)
1 /* Defines the basic SNES object */
2 #include <../src/snes/impls/fas/fasimpls.h>
3 
4 /*MC
5 Full Approximation Scheme nonlinear multigrid solver.
6 
7 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
8 
9 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
10 M*/
11 
12 extern PetscErrorCode SNESDestroy_FAS(SNES snes);
13 extern PetscErrorCode SNESSetUp_FAS(SNES snes);
14 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
15 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
16 extern PetscErrorCode SNESSolve_FAS(SNES snes);
17 extern PetscErrorCode SNESReset_FAS(SNES snes);
18 
19 EXTERN_C_BEGIN
20 
21 #undef __FUNCT__
22 #define __FUNCT__ "SNESCreate_FAS"
23 PetscErrorCode SNESCreate_FAS(SNES snes)
24 {
25   SNES_FAS * fas;
26   PetscErrorCode ierr;
27 
28   PetscFunctionBegin;
29   snes->ops->destroy        = SNESDestroy_FAS;
30   snes->ops->setup          = SNESSetUp_FAS;
31   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
32   snes->ops->view           = SNESView_FAS;
33   snes->ops->solve          = SNESSolve_FAS;
34   snes->ops->reset          = SNESReset_FAS;
35 
36   snes->usesksp             = PETSC_FALSE;
37   snes->usespc              = PETSC_FALSE;
38 
39   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
40   snes->data                = (void*) fas;
41   fas->level                = 0;
42   fas->levels               = 1;
43   fas->n_cycles             = 1;
44   fas->max_up_it            = 1;
45   fas->max_down_it          = 1;
46   fas->upsmooth             = PETSC_NULL;
47   fas->downsmooth           = PETSC_NULL;
48   fas->next                 = PETSC_NULL;
49   fas->interpolate          = PETSC_NULL;
50   fas->restrct              = PETSC_NULL;
51   fas->inject               = PETSC_NULL;
52 
53   PetscFunctionReturn(0);
54 }
55 EXTERN_C_END
56 
57 #undef __FUNCT__
58 #define __FUNCT__ "SNESFASGetLevels"
59 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
60   SNES_FAS * fas = (SNES_FAS *)snes->data;
61   PetscFunctionBegin;
62   *levels = fas->levels;
63   PetscFunctionReturn(0);
64 }
65 
66 #undef __FUNCT__
67 #define __FUNCT__ "SNESFASSetCycles"
68 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
69   SNES_FAS * fas = (SNES_FAS *)snes->data;
70   PetscErrorCode ierr;
71   PetscFunctionBegin;
72   fas->n_cycles = cycles;
73   if (fas->next) {
74     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
75   }
76   PetscFunctionReturn(0);
77 }
78 
79 
80 #undef __FUNCT__
81 #define __FUNCT__ "SNESFASSetCyclesOnLevel"
82 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, PetscInt cycles) {
83   SNES_FAS * fas =  (SNES_FAS *)snes->data;
84   PetscInt top_level = fas->level,i;
85 
86   PetscFunctionBegin;
87   if (level > top_level)
88     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
89   /* get to the correct level */
90   for (i = fas->level; i > level; i--) {
91     fas = (SNES_FAS *)fas->next->data;
92   }
93   if (fas->level != level)
94     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
95   fas->n_cycles = cycles;
96   PetscFunctionReturn(0);
97 }
98 
99 
100 #undef __FUNCT__
101 #define __FUNCT__ "SNESFASGetSNES"
102 PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
103   SNES_FAS * fas = (SNES_FAS *)snes->data;
104   PetscInt levels = fas->level;
105   PetscInt i;
106   PetscFunctionBegin;
107   *lsnes = snes;
108   if (fas->level < level) {
109     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
110   }
111   if (level > levels - 1) {
112     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
113   }
114   for (i = fas->level; i > level; i--) {
115     *lsnes = fas->next;
116     fas = (SNES_FAS *)(*lsnes)->data;
117   }
118   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
119   PetscFunctionReturn(0);
120 }
121 
122 #undef __FUNCT__
123 #define __FUNCT__ "SNESFASSetLevels"
124 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
125   PetscErrorCode ierr;
126   PetscInt i;
127   SNES_FAS * fas = (SNES_FAS *)snes->data;
128   MPI_Comm comm;
129   PetscFunctionBegin;
130   comm = ((PetscObject)snes)->comm;
131   if (levels == fas->levels) {
132     if (!comms)
133       PetscFunctionReturn(0);
134   }
135   /* user has changed the number of levels; reset */
136   ierr = SNESReset(snes);CHKERRQ(ierr);
137   /* destroy any coarser levels if necessary */
138   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
139   fas->next = PETSC_NULL;
140   /* setup the finest level */
141   for (i = levels - 1; i >= 0; i--) {
142     if (comms) comm = comms[i];
143     fas->level = i;
144     fas->levels = levels;
145     fas->next = PETSC_NULL;
146     if (i > 0) {
147       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
148       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
149       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
150       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
151       fas = (SNES_FAS *)fas->next->data;
152     }
153   }
154   PetscFunctionReturn(0);
155 }
156 
157 #undef __FUNCT__
158 #define __FUNCT__ "SNESFASSetInterpolation"
159 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
160   SNES_FAS * fas =  (SNES_FAS *)snes->data;
161   PetscInt top_level = fas->level,i;
162 
163   PetscFunctionBegin;
164   if (level > top_level)
165     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
166   /* get to the correct level */
167   for (i = fas->level; i > level; i--) {
168     fas = (SNES_FAS *)fas->next->data;
169   }
170   if (fas->level != level)
171     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
172   fas->interpolate = mat;
173   PetscFunctionReturn(0);
174 }
175 
176 #undef __FUNCT__
177 #define __FUNCT__ "SNESFASSetRestriction"
178 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
179   SNES_FAS * fas =  (SNES_FAS *)snes->data;
180   PetscInt top_level = fas->level,i;
181 
182   PetscFunctionBegin;
183   if (level > top_level)
184     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
185   /* get to the correct level */
186   for (i = fas->level; i > level; i--) {
187     fas = (SNES_FAS *)fas->next->data;
188   }
189   if (fas->level != level)
190     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
191   fas->restrct = mat;
192   PetscFunctionReturn(0);
193 }
194 
195 
196 #undef __FUNCT__
197 #define __FUNCT__ "SNESFASSetRScale"
198 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
199   SNES_FAS * fas =  (SNES_FAS *)snes->data;
200   PetscInt top_level = fas->level,i;
201 
202   PetscFunctionBegin;
203   if (level > top_level)
204     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
205   /* get to the correct level */
206   for (i = fas->level; i > level; i--) {
207     fas = (SNES_FAS *)fas->next->data;
208   }
209   if (fas->level != level)
210     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
211   fas->rscale = rscale;
212   PetscFunctionReturn(0);
213 }
214 
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "SNESFASSetInjection"
218 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
219   SNES_FAS * fas =  (SNES_FAS *)snes->data;
220   PetscInt top_level = fas->level,i;
221 
222   PetscFunctionBegin;
223   if (level > top_level)
224     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
225   /* get to the correct level */
226   for (i = fas->level; i > level; i--) {
227     fas = (SNES_FAS *)fas->next->data;
228   }
229   if (fas->level != level)
230     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
231   fas->inject = mat;
232   PetscFunctionReturn(0);
233 }
234 
235 #undef __FUNCT__
236 #define __FUNCT__ "SNESReset_FAS"
237 PetscErrorCode SNESReset_FAS(SNES snes)
238 {
239   PetscErrorCode ierr = 0;
240   SNES_FAS * fas = (SNES_FAS *)snes->data;
241 
242   PetscFunctionBegin;
243   /* destroy local data created in SNESSetup_FAS */
244 #if 0
245   /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */
246 #endif
247   if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr);
248 #if 0
249 #endif
250   PetscFunctionReturn(0);
251 }
252 
253 #undef __FUNCT__
254 #define __FUNCT__ "SNESDestroy_FAS"
255 PetscErrorCode SNESDestroy_FAS(SNES snes)
256 {
257   SNES_FAS * fas = (SNES_FAS *)snes->data;
258   PetscErrorCode ierr;
259 
260   PetscFunctionBegin;
261   /* recursively resets and then destroys */
262   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
263   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
264   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
265   if (fas->inject) {
266     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
267   }
268   if (fas->interpolate == fas->restrct) {
269     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
270     fas->restrct = PETSC_NULL;
271   } else {
272     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
273     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
274   }
275   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
276   if (snes->work)        ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);
277   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
278   ierr = PetscFree(fas);CHKERRQ(ierr);
279 
280   PetscFunctionReturn(0);
281 }
282 
283 #undef __FUNCT__
284 #define __FUNCT__ "SNESSetUp_FAS"
285 PetscErrorCode SNESSetUp_FAS(SNES snes)
286 {
287   SNES_FAS       *fas = (SNES_FAS *) snes->data,*tmp;
288   PetscErrorCode ierr;
289   VecScatter     injscatter;
290 
291   PetscFunctionBegin;
292   /* should call the SNESSetFromOptions() only when approriate */
293   tmp = fas;
294   while (tmp) {
295     if (tmp->upsmooth) {ierr = SNESSetFromOptions(tmp->upsmooth);CHKERRQ(ierr);}
296     if (tmp->downsmooth) {ierr = SNESSetFromOptions(tmp->downsmooth);CHKERRQ(ierr);}
297     tmp = tmp->next ? (SNES_FAS*) tmp->next->data : 0;
298   }
299 
300   if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */
301   /* gets the solver ready for solution */
302   if (fas->next) {
303     if (snes->ops->computegs) {
304       ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
305     }
306   }
307   if (snes->dm) {
308     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
309     if (fas->next) {
310       /* for now -- assume the DM and the evaluation functions have been set externally */
311       if (!fas->next->dm) {
312         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
313         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
314       }
315       /* set the interpolation and restriction from the DM */
316       if (!fas->interpolate) {
317         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
318         fas->restrct = fas->interpolate;
319       }
320       /* set the injection from the DM */
321       if (!fas->inject) {
322         ierr = DMGetInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
323         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
324         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
325       }
326     }
327     /* set the DMs of the pre and post-smoothers here */
328     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
329     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
330   }
331   if (fas->next) {
332     /* gotta set up the solution vector for this to work */
333     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
334     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
335   }
336   /* got to set them all up at once */
337   PetscFunctionReturn(0);
338 }
339 
340 #undef __FUNCT__
341 #define __FUNCT__ "SNESSetFromOptions_FAS"
342 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
343 {
344   SNES_FAS   *fas = (SNES_FAS *) snes->data;
345   PetscInt levels = 1;
346   PetscBool flg, monflg;
347   PetscErrorCode ierr;
348   const char * def_smooth = SNESNRICHARDSON;
349   char pre_type[256];
350   char post_type[256];
351   char                    monfilename[PETSC_MAX_PATH_LEN];
352 
353   PetscFunctionBegin;
354   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
355 
356   /* number of levels -- only process on the finest level */
357   if (fas->levels - 1 == fas->level) {
358     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
359     if (!flg && snes->dm) {
360       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
361       levels++;
362     }
363     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
364   }
365 
366   /* type of pre and/or post smoothers -- set both at once */
367   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
368   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
369   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
370   if (flg) {
371     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
372   } else {
373     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
374     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
375   }
376 
377   /* options for the number of preconditioning cycles and cycle type */
378   ierr = PetscOptionsInt("-snes_fas_up_it","Number of upsmoother iterations","PCMGSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
379   ierr = PetscOptionsInt("-snes_fas_down_it","Number of downsmoother iterations","PCMGSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
380   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","PCMGSetNumberSmoothUp",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
381 
382   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
383 
384   /* other options for the coarsest level */
385   if (fas->level == 0) {
386     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
387     if (flg) {
388       ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
389     } else {
390       ierr = PetscOptionsList("-snes_fas_coarse_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
391       ierr = PetscOptionsList("-snes_fas_coarse_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
392     }
393   }
394 
395   ierr = PetscOptionsTail();CHKERRQ(ierr);
396   /* setup from the determined types if the smoothers don't exist */
397   if (!fas->upsmooth) {
398     const char     *prefix;
399     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
400     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
401     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
402     if (fas->level || (fas->levels == 1)) {
403       ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_");CHKERRQ(ierr);
404     } else {
405       ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_coarse_");CHKERRQ(ierr);
406     }
407     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
408     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
409     if (snes->ops->computefunction) {
410       ierr = SNESSetFunction(fas->upsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
411     }
412   }
413   if (fas->upsmooth) {
414     ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
415   }
416 
417   if (!fas->downsmooth && fas->level != 0) {
418     const char     *prefix;
419     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
420     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
421     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
422     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_");CHKERRQ(ierr);
423     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
424     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
425     if (snes->ops->computefunction) {
426       ierr = SNESSetFunction(fas->downsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
427     }
428   }
429   if (fas->downsmooth) {
430     ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
431   }
432 
433   if (monflg) {
434     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
435 
436     /* set the monitors for the upsmoother and downsmoother */
437     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
438     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
439   }
440 
441   /* recursive option setting for the smoothers */
442   if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);}
443   PetscFunctionReturn(0);
444 }
445 
446 #undef __FUNCT__
447 #define __FUNCT__ "SNESView_FAS"
448 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
449 {
450   SNES_FAS   *fas = (SNES_FAS *) snes->data;
451   PetscBool      iascii;
452   PetscErrorCode ierr;
453   PetscInt levels = fas->levels;
454   PetscInt i;
455 
456   PetscFunctionBegin;
457   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
458   if (iascii) {
459     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
460     ierr = PetscViewerASCIIPushTab(viewer);
461     for (i = levels - 1; i >= 0; i--) {
462       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
463       if (fas->upsmooth) {
464         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
465         ierr = PetscViewerASCIIPushTab(viewer);
466         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
467         ierr = PetscViewerASCIIPopTab(viewer);
468       } else {
469         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
470       }
471       if (fas->downsmooth) {
472         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
473         ierr = PetscViewerASCIIPushTab(viewer);
474         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
475         ierr = PetscViewerASCIIPopTab(viewer);
476       } else {
477         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
478       }
479       if (fas->next) fas = (SNES_FAS *)fas->next->data;
480     }
481     ierr = PetscViewerASCIIPopTab(viewer);
482   } else {
483     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
484   }
485   PetscFunctionReturn(0);
486 }
487 
488 /*
489 
490 Defines the FAS cycle as:
491 
492 fine problem: F(x) = 0
493 coarse problem: F^c(x) = b^c
494 
495 b^c = F^c(I^c_fx^f - I^c_fF(x))
496 
497 correction:
498 
499 x = x + I(x^c - Rx)
500 
501  */
502 
503 #undef __FUNCT__
504 #define __FUNCT__ "FASCycle_Private"
505 PetscErrorCode FASCycle_Private(SNES snes, Vec X) {
506 
507   PetscErrorCode ierr;
508   Vec X_c, Xo_c, F_c, B_c,F,B;
509   SNES_FAS * fas = (SNES_FAS *)snes->data;
510   PetscInt i, k;
511   PetscScalar fnorm;
512 
513   PetscFunctionBegin;
514   F = snes->vec_func;
515   B = snes->vec_rhs;
516   /* pre-smooth -- just update using the pre-smoother */
517   if (snes->ops->computegs) {
518     if (fas->monitor) {
519       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
520       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
521       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
522       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d Gauss-Seidel Function norm %e\n", 0, fnorm);CHKERRQ(ierr);
523       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
524     }
525     for (k = 0; k < fas->max_up_it; k++) {
526       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
527       if (snes->monitor) {
528         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
529         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
530         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
531         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d Gauss-Seidel Function norm %e\n", k+1, fnorm);CHKERRQ(ierr);
532         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
533       }
534     }
535   } else if (fas->upsmooth) {
536     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
537   } else if (snes->pc) {
538     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
539   }
540   if (fas->next) {
541     for (i = 0; i < fas->n_cycles; i++) {
542       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
543 
544       X_c  = fas->next->vec_sol;
545       Xo_c = fas->next->work[0];
546       F_c  = fas->next->vec_func;
547       B_c  = fas->next->work[1];
548 
549       /* inject the solution */
550       if (fas->inject) {
551         ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
552       } else {
553         ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
554         ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
555       }
556       ierr = VecScale(F, -1.0);CHKERRQ(ierr);
557 
558       /* restrict the defect */
559       ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
560 
561       /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
562       fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
563       ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
564       ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
565 
566       /* set initial guess of the coarse problem to the projected fine solution */
567       ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
568 
569       /* recurse to the next level */
570       fas->next->vec_rhs = B_c;
571       ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr);
572       fas->next->vec_rhs = PETSC_NULL;
573 
574       /* correct as x <- x + I(x^c - Rx)*/
575       ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
576       ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr);
577     }
578   }
579     /* down-smooth -- just update using the down-smoother */
580   if (fas->level != 0) {
581     if (snes->ops->computegs) {
582       if (fas->monitor) {
583         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
584         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
585         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
586         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d Gauss-Seidel Function norm %e\n", 0, fnorm);CHKERRQ(ierr);
587         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
588       }
589       for (k = 0; k < fas->max_down_it; k++) {
590         ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
591         if (fas->monitor) {
592           ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
593           ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
594           ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
595           ierr = PetscViewerASCIIPrintf(fas->monitor, "%d Gauss-Seidel Function norm %e\n", k+1, fnorm);CHKERRQ(ierr);
596           ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
597         }
598       }
599     }else if (fas->downsmooth) {
600       ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
601     } else if (snes->pc) {
602       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
603     }
604   }
605   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
606 
607   PetscFunctionReturn(0);
608 }
609 
610 #undef __FUNCT__
611 #define __FUNCT__ "SNESSolve_FAS"
612 
613 PetscErrorCode SNESSolve_FAS(SNES snes)
614 {
615   PetscErrorCode ierr;
616   PetscInt i, maxits;
617   Vec X, B,F;
618   PetscReal fnorm;
619   PetscFunctionBegin;
620   maxits = snes->max_its;            /* maximum number of iterations */
621   snes->reason = SNES_CONVERGED_ITERATING;
622   X = snes->vec_sol;
623   B = snes->vec_rhs;
624   F = snes->vec_func;
625 
626   /*norm setup */
627   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
628   snes->iter = 0;
629   snes->norm = 0.;
630   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
631   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
632   if (snes->domainerror) {
633     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
634     PetscFunctionReturn(0);
635   }
636   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
637   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
638   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
639   snes->norm = fnorm;
640   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
641   SNESLogConvHistory(snes,fnorm,0);
642   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
643 
644   /* set parameter for default relative tolerance convergence test */
645   snes->ttol = fnorm*snes->rtol;
646   /* test convergence */
647   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
648   if (snes->reason) PetscFunctionReturn(0);
649   for (i = 0; i < maxits; i++) {
650     /* Call general purpose update function */
651 
652     if (snes->ops->update) {
653       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
654     }
655     ierr = FASCycle_Private(snes, X);CHKERRQ(ierr);
656     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
657     /* Monitor convergence */
658     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
659     snes->iter = i+1;
660     snes->norm = fnorm;
661     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
662     SNESLogConvHistory(snes,snes->norm,0);
663     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
664     /* Test for convergence */
665     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
666     if (snes->reason) break;
667   }
668   if (i == maxits) {
669     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
670     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
671   }
672   PetscFunctionReturn(0);
673 }
674