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