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