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