xref: /petsc/src/snes/impls/fas/fas.c (revision ab8d36c96b69fea5736263697ddaa972e87bcd95)
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","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
267   }
268 
269   ierr = PetscOptionsTail();CHKERRQ(ierr);
270   /* setup from the determined types if there is no pointwise procedure or smoother defined */
271 
272   if (monflg) {
273     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
274     /* set the monitors for the upsmoother and downsmoother */
275     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
276     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
277     if (fas->smoothu) ierr = SNESMonitorSet(fas->smoothu,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
278     if (fas->smoothd) ierr = SNESMonitorSet(fas->smoothd,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
279   } else {
280     /* unset the monitors on the coarse levels */
281     if (!isFine) {
282       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
283     }
284   }
285 
286   /* set up the default line search for coarse grid corrections */
287   if (fas->fastype == SNES_FAS_ADDITIVE) {
288     if (!snes->linesearch) {
289       ierr = SNESGetSNESLineSearch(snes, &linesearch);CHKERRQ(ierr);
290       ierr = SNESLineSearchSetType(linesearch, SNES_LINESEARCH_L2);CHKERRQ(ierr);
291     }
292   }
293 
294   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
295   /* recursive option setting for the smoothers */
296   if (next) {ierr = SNESSetFromOptions(next);CHKERRQ(ierr);}
297   PetscFunctionReturn(0);
298 }
299 
300 #undef __FUNCT__
301 #define __FUNCT__ "SNESView_FAS"
302 PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
303 {
304   SNES_FAS       *fas = (SNES_FAS *) snes->data;
305   PetscBool      isFine, iascii;
306   PetscInt       i;
307   PetscErrorCode ierr;
308   SNES           smoothu, smoothd, levelsnes;
309 
310   PetscFunctionBegin;
311   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
312   if (isFine) {
313     ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
314     if (iascii) {
315       ierr = PetscViewerASCIIPrintf(viewer, "FAS: type is %s, levels=%D, cycles=%D\n",  SNESFASTypes[fas->fastype], fas->levels, fas->n_cycles);CHKERRQ(ierr);
316       if (fas->galerkin) {
317         ierr = PetscViewerASCIIPrintf(viewer,"    Using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
318       } else {
319         ierr = PetscViewerASCIIPrintf(viewer,"    Not using Galerkin computed coarse grid function evaluation\n");CHKERRQ(ierr);
320       }
321       for (i=0; i<fas->levels; i++) {
322         ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
323         ierr = SNESFASCycleGetSmootherUp(levelsnes, &smoothu);CHKERRQ(ierr);
324         ierr = SNESFASCycleGetSmootherDown(levelsnes, &smoothd);CHKERRQ(ierr);
325         if (!i) {
326           ierr = PetscViewerASCIIPrintf(viewer,"Coarse grid solver -- level %D -------------------------------\n",i);CHKERRQ(ierr);
327         } else {
328           ierr = PetscViewerASCIIPrintf(viewer,"Down solver (pre-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
329         }
330         ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
331         ierr = SNESView(smoothd,viewer);CHKERRQ(ierr);
332         ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
333         if (i && (smoothd == smoothu)) {
334           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) same as down solver (pre-smoother)\n");CHKERRQ(ierr);
335         } else if (i) {
336           ierr = PetscViewerASCIIPrintf(viewer,"Up solver (post-smoother) on level %D -------------------------------\n",i);CHKERRQ(ierr);
337           ierr = PetscViewerASCIIPushTab(viewer);CHKERRQ(ierr);
338           ierr = SNESView(smoothu,viewer);CHKERRQ(ierr);
339           ierr = PetscViewerASCIIPopTab(viewer);CHKERRQ(ierr);
340         }
341       }
342     } else {
343       SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
344     }
345   }
346   PetscFunctionReturn(0);
347 }
348 
349 #undef __FUNCT__
350 #define __FUNCT__ "FASDownSmooth"
351 /*
352 Defines the action of the downsmoother
353  */
354 PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
355   PetscErrorCode      ierr = 0;
356   SNESConvergedReason reason;
357   Vec                 FPC;
358   SNES                smoothd;
359   PetscFunctionBegin;
360   ierr = SNESFASCycleGetSmootherDown(snes, &smoothd);CHKERRQ(ierr);
361   ierr = SNESSolve(smoothd, B, X);CHKERRQ(ierr);
362   /* check convergence reason for the smoother */
363   ierr = SNESGetConvergedReason(smoothd,&reason);CHKERRQ(ierr);
364   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
365     snes->reason = SNES_DIVERGED_INNER;
366     PetscFunctionReturn(0);
367   }
368   ierr = SNESGetFunction(smoothd, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
369   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
370   PetscFunctionReturn(0);
371 }
372 
373 
374 #undef __FUNCT__
375 #define __FUNCT__ "FASUpSmooth"
376 /*
377 Defines the action of the upsmoother
378  */
379 PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
380   PetscErrorCode      ierr = 0;
381   SNESConvergedReason reason;
382   Vec                 FPC;
383   SNES                smoothu;
384   PetscFunctionBegin;
385 
386   ierr = SNESFASCycleGetSmootherUp(snes, &smoothu);CHKERRQ(ierr);
387   ierr = SNESSolve(smoothu, B, X);CHKERRQ(ierr);
388   /* check convergence reason for the smoother */
389   ierr = SNESGetConvergedReason(smoothu,&reason);CHKERRQ(ierr);
390   if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
391     snes->reason = SNES_DIVERGED_INNER;
392     PetscFunctionReturn(0);
393   }
394   ierr = SNESGetFunction(smoothu, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
395   ierr = VecCopy(FPC, F);CHKERRQ(ierr);
396   PetscFunctionReturn(0);
397 }
398 
399 #undef __FUNCT__
400 #define __FUNCT__ "SNESFASCreateCoarseVec"
401 /*@
402    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
403 
404    Collective
405 
406    Input Arguments:
407 .  snes - SNESFAS
408 
409    Output Arguments:
410 .  Xcoarse - vector on level one coarser than snes
411 
412    Level: developer
413 
414 .seealso: SNESFASSetRestriction(), SNESFASRestrict()
415 @*/
416 PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
417 {
418   PetscErrorCode ierr;
419   SNES_FAS       *fas = (SNES_FAS*)snes->data;
420 
421   PetscFunctionBegin;
422   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
423   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
424   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
425   PetscFunctionReturn(0);
426 }
427 
428 #undef __FUNCT__
429 #define __FUNCT__ "SNESFASRestrict"
430 /*@
431    SNESFASRestrict - restrict a Vec to the next coarser level
432 
433    Collective
434 
435    Input Arguments:
436 +  fine - SNES from which to restrict
437 -  Xfine - vector to restrict
438 
439    Output Arguments:
440 .  Xcoarse - result of restriction
441 
442    Level: developer
443 
444 .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
445 @*/
446 PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
447 {
448   PetscErrorCode ierr;
449   SNES_FAS       *fas = (SNES_FAS*)fine->data;
450 
451   PetscFunctionBegin;
452   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
453   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
454   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
455   if (fas->inject) {
456     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
457   } else {
458     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
459     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
460   }
461   PetscFunctionReturn(0);
462 }
463 
464 #undef __FUNCT__
465 #define __FUNCT__ "FASCoarseCorrection"
466 /*
467 
468 Performs the FAS coarse correction as:
469 
470 fine problem: F(x) = 0
471 coarse problem: F^c(x) = b^c
472 
473 b^c = F^c(I^c_fx^f - I^c_fF(x))
474 
475  */
476 PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
477   PetscErrorCode      ierr;
478   Vec                 X_c, Xo_c, F_c, B_c;
479   SNESConvergedReason reason;
480   SNES                next;
481   Mat                 restrct, interpolate;
482   PetscFunctionBegin;
483   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
484   if (next) {
485     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
486     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
487 
488     X_c  = next->vec_sol;
489     Xo_c = next->work[0];
490     F_c  = next->vec_func;
491     B_c  = next->vec_rhs;
492 
493     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
494     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
495 
496     /* restrict the defect */
497     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
498 
499     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
500     next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
501     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
502 
503     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
504 
505     /* set initial guess of the coarse problem to the projected fine solution */
506     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
507 
508     /* recurse to the next level */
509     next->vec_rhs = B_c;
510     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
511     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
512     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
513       snes->reason = SNES_DIVERGED_INNER;
514       PetscFunctionReturn(0);
515     }
516 
517     /* correct as x <- x + I(x^c - Rx)*/
518     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
519     ierr = MatInterpolateAdd(interpolate, X_c, X, X_new);CHKERRQ(ierr);
520   }
521   PetscFunctionReturn(0);
522 }
523 
524 #undef __FUNCT__
525 #define __FUNCT__ "FASCycle_Additive"
526 /*
527 
528 The additive cycle looks like:
529 
530 xhat = x
531 xhat = dS(x, b)
532 x = coarsecorrection(xhat, b_d)
533 x = x + nu*(xhat - x);
534 (optional) x = uS(x, b)
535 
536 With the coarse RHS (defect correction) as below.
537 
538  */
539 PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
540   Vec                 F, B, Xhat;
541   Vec                 X_c, Xo_c, F_c, B_c;
542   PetscErrorCode      ierr;
543   SNESConvergedReason reason;
544   PetscReal           xnorm, fnorm, ynorm;
545   PetscBool           lssuccess;
546   SNES                next;
547   Mat                 restrct, interpolate;
548   PetscFunctionBegin;
549   ierr = SNESFASCycleGetCorrection(snes, &next);CHKERRQ(ierr);
550   F = snes->vec_func;
551   B = snes->vec_rhs;
552   Xhat = snes->work[3];
553   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
554   /* recurse first */
555   if (next) {
556     ierr = SNESFASCycleGetRestriction(snes, &restrct);CHKERRQ(ierr);
557     ierr = SNESFASCycleGetInterpolation(snes, &interpolate);CHKERRQ(ierr);
558     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
559 
560     X_c  = next->vec_sol;
561     Xo_c = next->work[0];
562     F_c  = next->vec_func;
563     B_c  = next->vec_rhs;
564 
565     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
566     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
567 
568     /* restrict the defect */
569     ierr = MatRestrict(restrct, F, B_c);CHKERRQ(ierr);
570 
571     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
572     next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
573     ierr = SNESComputeFunction(next, Xo_c, F_c);CHKERRQ(ierr);
574 
575     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
576 
577     /* set initial guess of the coarse problem to the projected fine solution */
578     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
579 
580     /* recurse */
581     next->vec_rhs = B_c;
582     ierr = SNESSolve(next, B_c, X_c);CHKERRQ(ierr);
583 
584     /* smooth on this level */
585     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
586 
587     ierr = SNESGetConvergedReason(next,&reason);CHKERRQ(ierr);
588     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
589       snes->reason = SNES_DIVERGED_INNER;
590       PetscFunctionReturn(0);
591     }
592 
593     /* correct as x <- x + I(x^c - Rx)*/
594     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
595     ierr = MatInterpolate(interpolate, X_c, Xhat);CHKERRQ(ierr);
596 
597     /* additive correction of the coarse direction*/
598     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
599     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
600     ierr = SNESLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
601     ierr = SNESLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
602     if (!lssuccess) {
603       if (++snes->numFailures >= snes->maxFailures) {
604         snes->reason = SNES_DIVERGED_LINE_SEARCH;
605         PetscFunctionReturn(0);
606       }
607     }
608     ierr = SNESLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
609   } else {
610     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
611   }
612   PetscFunctionReturn(0);
613 }
614 
615 #undef __FUNCT__
616 #define __FUNCT__ "FASCycle_Multiplicative"
617 /*
618 
619 Defines the FAS cycle as:
620 
621 fine problem: F(x) = 0
622 coarse problem: F^c(x) = b^c
623 
624 b^c = F^c(I^c_fx^f - I^c_fF(x))
625 
626 correction:
627 
628 x = x + I(x^c - Rx)
629 
630  */
631 PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
632 
633   PetscErrorCode      ierr;
634   Vec                 F,B;
635   SNES_FAS            *fas = (SNES_FAS *)snes->data;
636 
637   PetscFunctionBegin;
638   F = snes->vec_func;
639   B = snes->vec_rhs;
640   /* pre-smooth -- just update using the pre-smoother */
641   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
642 
643   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
644 
645   if (fas->level != 0) {
646     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
647   }
648   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
649 
650   PetscFunctionReturn(0);
651 }
652 
653 #undef __FUNCT__
654 #define __FUNCT__ "SNESSolve_FAS"
655 
656 PetscErrorCode SNESSolve_FAS(SNES snes)
657 {
658   PetscErrorCode ierr;
659   PetscInt       i, maxits;
660   Vec            X, F;
661   PetscReal      fnorm;
662   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
663   DM             dm;
664 
665   PetscFunctionBegin;
666   maxits = snes->max_its;            /* maximum number of iterations */
667   snes->reason = SNES_CONVERGED_ITERATING;
668   X = snes->vec_sol;
669   F = snes->vec_func;
670 
671   /*norm setup */
672   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
673   snes->iter = 0;
674   snes->norm = 0.;
675   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
676   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
677   if (snes->domainerror) {
678     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
679     PetscFunctionReturn(0);
680   }
681   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
682   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
683   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
684   snes->norm = fnorm;
685   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
686   SNESLogConvHistory(snes,fnorm,0);
687   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
688 
689   /* set parameter for default relative tolerance convergence test */
690   snes->ttol = fnorm*snes->rtol;
691   /* test convergence */
692   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
693   if (snes->reason) PetscFunctionReturn(0);
694 
695   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
696   for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
697     DM dmcoarse;
698     ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
699     ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
700     dm = dmcoarse;
701   }
702 
703   for (i = 0; i < maxits; i++) {
704     /* Call general purpose update function */
705 
706     if (snes->ops->update) {
707       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
708     }
709     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
710       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
711     } else {
712       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
713     }
714 
715     /* check for FAS cycle divergence */
716     if (snes->reason != SNES_CONVERGED_ITERATING) {
717       PetscFunctionReturn(0);
718     }
719     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
720     /* Monitor convergence */
721     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
722     snes->iter = i+1;
723     snes->norm = fnorm;
724     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
725     SNESLogConvHistory(snes,snes->norm,0);
726     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
727     /* Test for convergence */
728     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
729     if (snes->reason) break;
730   }
731   if (i == maxits) {
732     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
733     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
734   }
735   PetscFunctionReturn(0);
736 }
737