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