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