xref: /petsc/src/snes/impls/fas/fas.c (revision c8c899ca14d28eb40f2b95e79d22697c89496ce1)
1 /* Defines the basic SNES object */
2 #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3 
4 const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0};
5 
6 extern PetscErrorCode SNESDestroy_FAS(SNES snes);
7 extern PetscErrorCode SNESSetUp_FAS(SNES snes);
8 extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
9 extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
10 extern PetscErrorCode SNESSolve_FAS(SNES snes);
11 extern PetscErrorCode SNESReset_FAS(SNES snes);
12 extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *);
13 extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *);
14 
15 /*MC
16 
17 SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
18 
19 The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
20 
21 Level: advanced
22 
23 .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
24 M*/
25 
26 #undef __FUNCT__
27 #define __FUNCT__ "SNESCreate_FAS"
28 PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes)
29 {
30   SNES_FAS * fas;
31   PetscErrorCode ierr;
32 
33   PetscFunctionBegin;
34   snes->ops->destroy        = SNESDestroy_FAS;
35   snes->ops->setup          = SNESSetUp_FAS;
36   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
37   snes->ops->view           = SNESView_FAS;
38   snes->ops->solve          = SNESSolve_FAS;
39   snes->ops->reset          = SNESReset_FAS;
40 
41   snes->usesksp             = PETSC_FALSE;
42   snes->usespc              = PETSC_FALSE;
43 
44   if (!snes->tolerancesset) {
45     snes->max_funcs = 30000;
46     snes->max_its   = 10000;
47   }
48 
49   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
50   snes->data                  = (void*) fas;
51   fas->level                  = 0;
52   fas->levels                 = 1;
53   fas->n_cycles               = 1;
54   fas->max_up_it              = 1;
55   fas->max_down_it            = 1;
56   fas->smoothu                = PETSC_NULL;
57   fas->smoothd                = PETSC_NULL;
58   fas->next                   = PETSC_NULL;
59   fas->previous               = PETSC_NULL;
60   fas->fine                   = snes;
61   fas->interpolate            = PETSC_NULL;
62   fas->restrct                = PETSC_NULL;
63   fas->inject                 = PETSC_NULL;
64   fas->monitor                = PETSC_NULL;
65   fas->usedmfornumberoflevels = PETSC_FALSE;
66   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
67   PetscFunctionReturn(0);
68 }
69 
70 #undef __FUNCT__
71 #define __FUNCT__ "SNESReset_FAS"
72 PetscErrorCode SNESReset_FAS(SNES snes)
73 {
74   PetscErrorCode ierr = 0;
75   SNES_FAS * fas = (SNES_FAS *)snes->data;
76 
77   PetscFunctionBegin;
78   if (fas->smoothu != fas->smoothd) {
79     ierr = SNESDestroy(&fas->smoothu);CHKERRQ(ierr);
80     ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
81   } else {
82     ierr = SNESDestroy(&fas->smoothd);CHKERRQ(ierr);
83     fas->smoothu = PETSC_NULL;
84   }
85   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
86   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
87   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
88   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
89   if (fas->next)      ierr = SNESReset(fas->next);CHKERRQ(ierr);
90 
91   PetscFunctionReturn(0);
92 }
93 
94 #undef __FUNCT__
95 #define __FUNCT__ "SNESDestroy_FAS"
96 PetscErrorCode SNESDestroy_FAS(SNES snes)
97 {
98   SNES_FAS * fas = (SNES_FAS *)snes->data;
99   PetscErrorCode ierr = 0;
100 
101   PetscFunctionBegin;
102   /* recursively resets and then destroys */
103   ierr = SNESReset(snes);CHKERRQ(ierr);
104   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
105   ierr = PetscFree(fas);CHKERRQ(ierr);
106 
107   PetscFunctionReturn(0);
108 }
109 
110 #undef __FUNCT__
111 #define __FUNCT__ "SNESSetUp_FAS"
112 PetscErrorCode SNESSetUp_FAS(SNES snes)
113 {
114   SNES_FAS                *fas = (SNES_FAS *) snes->data;
115   PetscErrorCode          ierr;
116   VecScatter              injscatter;
117   PetscInt                dm_levels;
118   Vec                     vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
119   SNES                    next;
120   PetscBool               isFine;
121   PetscFunctionBegin;
122 
123   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
124   if (fas->usedmfornumberoflevels && isFine) {
125     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
126     dm_levels++;
127     if (dm_levels > fas->levels) {
128       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
129       vec_sol = snes->vec_sol;
130       vec_func = snes->vec_func;
131       vec_sol_update = snes->vec_sol_update;
132       vec_rhs = snes->vec_rhs;
133       snes->vec_sol = PETSC_NULL;
134       snes->vec_func = PETSC_NULL;
135       snes->vec_sol_update = PETSC_NULL;
136       snes->vec_rhs = PETSC_NULL;
137 
138       /* reset the number of levels */
139       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
140       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
141 
142       snes->vec_sol = vec_sol;
143       snes->vec_func = vec_func;
144       snes->vec_rhs = vec_rhs;
145       snes->vec_sol_update = vec_sol_update;
146     }
147   }
148   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
149   if (!isFine) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
150 
151   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
152     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
153   } else {
154     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
155   }
156 
157   /* set up the smoothers if they haven't already been set up */
158   if (!fas->smoothd) {
159     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
160   }
161 
162   if (snes->dm) {
163     /* set the smoother DMs properly */
164     if (fas->smoothu) ierr = SNESSetDM(fas->smoothu, snes->dm);CHKERRQ(ierr);
165     ierr = SNESSetDM(fas->smoothd, snes->dm);CHKERRQ(ierr);
166     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
167     if (next) {
168       /* for now -- assume the DM and the evaluation functions have been set externally */
169       if (!next->dm) {
170         ierr = DMCoarsen(snes->dm, ((PetscObject)next)->comm, &next->dm);CHKERRQ(ierr);
171         ierr = SNESSetDM(next, next->dm);CHKERRQ(ierr);
172       }
173       /* set the interpolation and restriction from the DM */
174       if (!fas->interpolate) {
175         ierr = DMCreateInterpolation(next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
176         if (!fas->restrct) {
177           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
178           fas->restrct = fas->interpolate;
179         }
180       }
181       /* set the injection from the DM */
182       if (!fas->inject) {
183         ierr = DMCreateInjection(next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
184         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
185         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
186       }
187     }
188   }
189   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
190 
191   if (fas->galerkin) {
192     if (next)
193       ierr = SNESSetFunction(next, PETSC_NULL, SNESFASGalerkinDefaultFunction, next);CHKERRQ(ierr);
194     if (fas->smoothd && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothd, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
195     if (fas->smoothu && fas->level != fas->levels - 1) ierr = SNESSetFunction(fas->smoothu, PETSC_NULL, SNESFASGalerkinDefaultFunction, snes);CHKERRQ(ierr);
196   }
197 
198   if (next) {
199     /* gotta set up the solution vector for this to work */
200     if (!next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_sol);CHKERRQ(ierr);}
201     if (!next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&next->vec_rhs);CHKERRQ(ierr);}
202     ierr = SNESSetUp(next);CHKERRQ(ierr);
203   }
204   /* setup FAS work vectors */
205   if (fas->galerkin) {
206     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
207     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
208   }
209   PetscFunctionReturn(0);
210 }
211 
212 #undef __FUNCT__
213 #define __FUNCT__ "SNESSetFromOptions_FAS"
214 PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
215 {
216   SNES_FAS       *fas = (SNES_FAS *) snes->data;
217   PetscInt       levels = 1;
218   PetscBool      flg, monflg, galerkinflg;
219   PetscErrorCode ierr;
220   char           monfilename[PETSC_MAX_PATH_LEN];
221   SNESFASType    fastype;
222   const char     *optionsprefix;
223   SNESLineSearch linesearch;
224   PetscInt       m;
225   SNES           next;
226   PetscBool      isFine;
227 
228   PetscFunctionBegin;
229   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
230   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
231 
232   /* number of levels -- only process most options on the finest level */
233   if (isFine) {
234     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
235     if (!flg && snes->dm) {
236       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
237       levels++;
238       fas->usedmfornumberoflevels = PETSC_TRUE;
239     }
240     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
241     fastype = fas->fastype;
242     ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
243     if (flg) {
244       ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
245     }
246 
247     ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
248     ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&m,&flg);CHKERRQ(ierr);
249     if (flg) {
250       ierr = SNESFASSetCycles(snes, m);CHKERRQ(ierr);
251     }
252 
253     ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFASSetGalerkin",fas->galerkin,&galerkinflg,&flg);CHKERRQ(ierr);
254     if (flg) {
255       ierr = SNESFASSetGalerkin(snes, galerkinflg);CHKERRQ(ierr);
256     }
257 
258     ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smoothing steps","SNESFASSetNumberSmoothUp",1,&m,&flg);CHKERRQ(ierr);
259     if (flg) {
260       ierr = SNESFASSetNumberSmoothUp(snes,m);CHKERRQ(ierr);
261     }
262     ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smoothing steps","SNESFASSetNumberSmoothDown",1,&m,&flg);CHKERRQ(ierr);
263     if (flg) {
264       ierr = SNESFASSetNumberSmoothDown(snes,m);CHKERRQ(ierr);
265     }
266     ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESFASSetMonitor","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
267     if (monflg) ierr = SNESFASSetMonitor(snes, PETSC_TRUE);CHKERRQ(ierr);
268   }
269 
270   ierr = PetscOptionsTail();CHKERRQ(ierr);
271   /* setup from the determined types if there is no pointwise procedure or smoother defined */
272 
273   /* set up the default line search for coarse grid corrections */
274   if (fas->fastype == SNES_FAS_ADDITIVE) {
275     if (!snes->linesearch) {
276       ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
277       ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr);
278     }
279   }
280 
281   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
282   /* recursive option setting for the smoothers */
283   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
284   PetscFunctionReturn(0);
285 }
286 
287 #undef __FUNCT__
288 #define __FUNCT__ "SNESView_FAS"
289 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
290 {
291   SNES_FAS       *fas = (SNES_FAS *) snes->data;
292   PetscBool      isFine, iascii;
293   PetscInt       i;
294   PetscErrorCode ierr;
295   SNES           smoothu, smoothd, levelsnes;
296 
297   PetscFunctionBegin;
298   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
299   if (isFine) {
300     ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
301     if (iascii) {
302       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
303       if (fas->galerkin) {
304         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
305       } else {
306         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
307       }
308       for (i=0; i<fas->levels; i++) {
309         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
310         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
311         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
312         if (!i) {
313           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
314         } else {
315           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
316         }
317         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
318         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
319         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
320         if (i && (smoothd == smoothu)) {
321           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
322         } else if (i) {
323           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
324           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
325           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
326           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
327         }
328       }
329     } else {
330       SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
331     }
332   }
333   PetscFunctionReturn(0);
334 }
335 
336 #undef __FUNCT__
337 #define __FUNCT__ "FASDownSmooth"
338 /*
339 Defines the action of the downsmoother
340  */
341 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
342   PetscErrorCode      ierr = 0;
343   SNESConvergedReason reason;
344   Vec                 FPC;
345   SNES                smoothd;
346   PetscFunctionBegin;
347   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
348   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
349   /* check convergence reason for the smoother */
350   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
351   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
352     snes->reason = SNES_DIVERGED_INNER;
353     PetscFunctionReturn(0);
354   }
355   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
356   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
357   PetscFunctionReturn(0);
358 }
359 
360 
361 #undef __FUNCT__
362 #define __FUNCT__ "FASUpSmooth"
363 /*
364 Defines the action of the upsmoother
365  */
366 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
367   PetscErrorCode      ierr = 0;
368   SNESConvergedReason reason;
369   Vec                 FPC;
370   SNES                smoothu;
371   PetscFunctionBegin;
372 
373   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
374   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
375   /* check convergence reason for the smoother */
376   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
377   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
378     snes->reason = SNES_DIVERGED_INNER;
379     PetscFunctionReturn(0);
380   }
381   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
382   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
383   PetscFunctionReturn(0);
384 }
385 
386 #undef __FUNCT__
387 #define __FUNCT__ "SNESFASCreateCoarseVec"
388 /*@
389    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
390 
391    Collective
392 
393    Input Arguments:
394 .  snes - SNESFAS
395 
396    Output Arguments:
397 .  Xcoarse - vector on level one coarser than snes
398 
399    Level: developer
400 
401 .seealso: SNESFASSetRestriction(), SNESFASRestrict()
402 @*/
403 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
404 {
405   PetscErrorCode ierr;
406   SNES_FAS       *fas = (SNES_FAS*)snes->data;
407 
408   PetscFunctionBegin;
409   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
410   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
411   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
412   PetscFunctionReturn(0);
413 }
414 
415 #undef __FUNCT__
416 #define __FUNCT__ "SNESFASRestrict"
417 /*@
418    SNESFASRestrict - restrict a Vec to the next coarser level
419 
420    Collective
421 
422    Input Arguments:
423 +  fine - SNES from which to restrict
424 -  Xfine - vector to restrict
425 
426    Output Arguments:
427 .  Xcoarse - result of restriction
428 
429    Level: developer
430 
431 .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
432 @*/
433 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
434 {
435   PetscErrorCode ierr;
436   SNES_FAS       *fas = (SNES_FAS*)fine->data;
437 
438   PetscFunctionBegin;
439   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
440   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
441   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
442   if (fas->inject) {
443     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
444   } else {
445     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
446     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
447   }
448   PetscFunctionReturn(0);
449 }
450 
451 #undef __FUNCT__
452 #define __FUNCT__ "FASCoarseCorrection"
453 /*
454 
455 Performs the FAS coarse correction as:
456 
457 fine problem: F(x) = 0
458 coarse problem: F^c(x) = b^c
459 
460 b^c = F^c(I^c_fx^f - I^c_fF(x))
461 
462  */
463 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
464   PetscErrorCode      ierr;
465   Vec                 X_c, Xo_c, F_c, B_c;
466   SNESConvergedReason reason;
467   SNES                next;
468   Mat                 restrct, interpolate;
469   PetscFunctionBegin;
470   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
471   if (next) {
472     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
473     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
474 
475     X_c  = next->vec_sol;
476     Xo_c = next->work[0];
477     F_c  = next->vec_func;
478     B_c  = next->vec_rhs;
479 
480     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
481     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
482 
483     /* restrict the defect */
484     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
485 
486     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
487     next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
488     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
489 
490     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
491 
492     /* set initial guess of the coarse problem to the projected fine solution */
493     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
494 
495     /* recurse to the next level */
496     next->vec_rhs = B_c;
497     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
498     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
499     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
500       snes->reason = SNES_DIVERGED_INNER;
501       PetscFunctionReturn(0);
502     }
503 
504     /* correct as x <- x + I(x^c - Rx)*/
505     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
506     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
507   }
508   PetscFunctionReturn(0);
509 }
510 
511 #undef __FUNCT__
512 #define __FUNCT__ "FASCycle_Additive"
513 /*
514 
515 The additive cycle looks like:
516 
517 xhat = x
518 xhat = dS(x, b)
519 x = coarsecorrection(xhat, b_d)
520 x = x + nu*(xhat - x);
521 (optional) x = uS(x, b)
522 
523 With the coarse RHS (defect correction) as below.
524 
525  */
526 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
527   Vec                 F, B, Xhat;
528   Vec                 X_c, Xo_c, F_c, B_c;
529   PetscErrorCode      ierr;
530   SNESConvergedReason reason;
531   PetscReal           xnorm, fnorm, ynorm;
532   PetscBool           lssuccess;
533   SNES                next;
534   Mat                 restrct, interpolate;
535   PetscFunctionBegin;
536   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
537   F = snes->vec_func;
538   B = snes->vec_rhs;
539   Xhat = snes->work[3];
540   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
541   /* recurse first */
542   if (next) {
543     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
544     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
545     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
546 
547     X_c  = next->vec_sol;
548     Xo_c = next->work[0];
549     F_c  = next->vec_func;
550     B_c  = next->vec_rhs;
551 
552     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
553     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
554 
555     /* restrict the defect */
556     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
557 
558     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
559     next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
560     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
561 
562     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
563 
564     /* set initial guess of the coarse problem to the projected fine solution */
565     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
566 
567     /* recurse */
568     next->vec_rhs = B_c;
569     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
570 
571     /* smooth on this level */
572     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
573 
574     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
575     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
576       snes->reason = SNES_DIVERGED_INNER;
577       PetscFunctionReturn(0);
578     }
579 
580     /* correct as x <- x + I(x^c - Rx)*/
581     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
582     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
583 
584     /* additive correction of the coarse direction*/
585     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
586     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
587     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
588     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
589     if (!lssuccess) {
590       if (++snes->numFailures >= snes->maxFailures) {
591         snes->reason = SNES_DIVERGED_LINE_SEARCH;
592         PetscFunctionReturn(0);
593       }
594     }
595     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
596   } else {
597     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
598   }
599   PetscFunctionReturn(0);
600 }
601 
602 #undef __FUNCT__
603 #define __FUNCT__ "FASCycle_Multiplicative"
604 /*
605 
606 Defines the FAS cycle as:
607 
608 fine problem: F(x) = 0
609 coarse problem: F^c(x) = b^c
610 
611 b^c = F^c(I^c_fx^f - I^c_fF(x))
612 
613 correction:
614 
615 x = x + I(x^c - Rx)
616 
617  */
618 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
619 
620   PetscErrorCode      ierr;
621   Vec                 F,B;
622   SNES_FAS            *fas = (SNES_FAS *)snes->data;
623 
624   PetscFunctionBegin;
625   F = snes->vec_func;
626   B = snes->vec_rhs;
627   /* pre-smooth -- just update using the pre-smoother */
628   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
629 
630   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
631 
632   if (fas->level != 0) {
633     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
634   }
635   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
636 
637   PetscFunctionReturn(0);
638 }
639 
640 #undef __FUNCT__
641 #define __FUNCT__ "SNESSolve_FAS"
642 
643 PetscErrorCode SNESSolve_FAS(SNES snes)
644 {
645   PetscErrorCode ierr;
646   PetscInt       i, maxits;
647   Vec            X, F;
648   PetscReal      fnorm;
649   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
650   DM             dm;
651 
652   PetscFunctionBegin;
653   maxits = snes->max_its;            /* maximum number of iterations */
654   snes->reason = SNES_CONVERGED_ITERATING;
655   X = snes->vec_sol;
656   F = snes->vec_func;
657 
658   /*norm setup */
659   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
660   snes->iter = 0;
661   snes->norm = 0.;
662   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
663   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
664   if (snes->domainerror) {
665     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
666     PetscFunctionReturn(0);
667   }
668   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
669   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
670   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
671   snes->norm = fnorm;
672   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
673   SNESLogConvHistory(snes,fnorm,0);
674   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
675 
676   /* set parameter for default relative tolerance convergence test */
677   snes->ttol = fnorm*snes->rtol;
678   /* test convergence */
679   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
680   if (snes->reason) PetscFunctionReturn(0);
681 
682   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
683   for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
684     DM dmcoarse;
685     ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
686     ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
687     dm = dmcoarse;
688   }
689 
690   for (i = 0; i < maxits; i++) {
691     /* Call general purpose update function */
692 
693     if (snes->ops->update) {
694       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
695     }
696     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
697       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
698     } else {
699       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
700     }
701 
702     /* check for FAS cycle divergence */
703     if (snes->reason != SNES_CONVERGED_ITERATING) {
704       PetscFunctionReturn(0);
705     }
706     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
707     /* Monitor convergence */
708     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
709     snes->iter = i+1;
710     snes->norm = fnorm;
711     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
712     SNESLogConvHistory(snes,snes->norm,0);
713     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
714     /* Test for convergence */
715     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
716     if (snes->reason) break;
717   }
718   if (i == maxits) {
719     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
720     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
721   }
722   PetscFunctionReturn(0);
723 }
724