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