xref: /petsc/src/snes/impls/fas/fas.c (revision fa9694d7f50f79d530566ef7a607d5b0aa6ae5fc)
1421d9b32SPeter Brune /* Defines the basic SNES object */
2421d9b32SPeter Brune #include <../src/snes/impls/fas/fasimpls.h>
3421d9b32SPeter Brune 
4421d9b32SPeter Brune /*MC
5421d9b32SPeter Brune Full Approximation Scheme nonlinear multigrid solver.
6421d9b32SPeter Brune 
7421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
8421d9b32SPeter Brune 
9421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
10421d9b32SPeter Brune M*/
11421d9b32SPeter Brune 
12421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
13421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
14421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
15421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
16421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
17421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
18421d9b32SPeter Brune 
19421d9b32SPeter Brune EXTERN_C_BEGIN
20421d9b32SPeter Brune 
21421d9b32SPeter Brune #undef __FUNCT__
22421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
23421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes)
24421d9b32SPeter Brune {
25421d9b32SPeter Brune   SNES_FAS * fas;
26421d9b32SPeter Brune   PetscErrorCode ierr;
27421d9b32SPeter Brune 
28421d9b32SPeter Brune   PetscFunctionBegin;
29421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
30421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
31421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
32421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
33421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
34421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
35421d9b32SPeter Brune 
36421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
37421d9b32SPeter Brune   snes->data = (void*) fas;
38421d9b32SPeter Brune   fas->level = 0;
39421d9b32SPeter Brune   fas->presmooth  = PETSC_NULL;
40421d9b32SPeter Brune   fas->postsmooth = PETSC_NULL;
41421d9b32SPeter Brune   fas->next = PETSC_NULL;
42421d9b32SPeter Brune   fas->interpolate = PETSC_NULL;
43421d9b32SPeter Brune   fas->restrct = PETSC_NULL;
44421d9b32SPeter Brune 
45421d9b32SPeter Brune   PetscFunctionReturn(0);
46421d9b32SPeter Brune }
47421d9b32SPeter Brune EXTERN_C_END
48421d9b32SPeter Brune 
49421d9b32SPeter Brune #undef __FUNCT__
50421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
51421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
52421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
53421d9b32SPeter Brune   PetscFunctionBegin;
54421d9b32SPeter Brune   *levels = fas->level;
55421d9b32SPeter Brune   PetscFunctionReturn(0);
56421d9b32SPeter Brune }
57421d9b32SPeter Brune 
58421d9b32SPeter Brune #undef __FUNCT__
59421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
60421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
61421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
62421d9b32SPeter Brune   PetscInt levels = fas->level;
63421d9b32SPeter Brune   PetscInt i;
64421d9b32SPeter Brune   PetscFunctionBegin;
65421d9b32SPeter Brune   *lsnes = snes;
66421d9b32SPeter Brune   if (fas->level < level) {
67421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
68421d9b32SPeter Brune   }
69421d9b32SPeter Brune   if (level > levels - 1) {
70421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
71421d9b32SPeter Brune   }
72421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
73421d9b32SPeter Brune     *lsnes = fas->next;
74421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
75421d9b32SPeter Brune   }
76421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
77421d9b32SPeter Brune   PetscFunctionReturn(0);
78421d9b32SPeter Brune }
79421d9b32SPeter Brune 
80421d9b32SPeter Brune #undef __FUNCT__
81421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
82421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
83421d9b32SPeter Brune   PetscErrorCode ierr;
84421d9b32SPeter Brune   PetscInt i;
85421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
86421d9b32SPeter Brune   MPI_Comm comm;
87421d9b32SPeter Brune 
88421d9b32SPeter Brune   PetscFunctionBegin;
89421d9b32SPeter Brune   comm = PETSC_COMM_WORLD;
90421d9b32SPeter Brune   /* user has changed the number of levels; reset */
91421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
92421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
93421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
94421d9b32SPeter Brune   /* setup the finest level */
95421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
96421d9b32SPeter Brune     if (comms) comm = comms[i];
97421d9b32SPeter Brune     fas->level = i;
98421d9b32SPeter Brune     fas->levels = levels;
99421d9b32SPeter Brune     if (i > 0) {
100421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
101421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
102421d9b32SPeter Brune       fas = (SNES_FAS *)fas->next->data;
103421d9b32SPeter Brune     }
104421d9b32SPeter Brune   }
105421d9b32SPeter Brune   PetscFunctionReturn(0);
106421d9b32SPeter Brune }
107421d9b32SPeter Brune 
108421d9b32SPeter Brune #undef __FUNCT__
109421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
110421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
111421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
112d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
113421d9b32SPeter Brune 
114421d9b32SPeter Brune   PetscFunctionBegin;
115421d9b32SPeter Brune   if (level > top_level)
116421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
117421d9b32SPeter Brune   /* get to the correct level */
118d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
119421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
120421d9b32SPeter Brune   }
121421d9b32SPeter Brune   if (fas->level != level)
122421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
123421d9b32SPeter Brune   fas->interpolate = mat;
124421d9b32SPeter Brune   PetscFunctionReturn(0);
125421d9b32SPeter Brune }
126421d9b32SPeter Brune 
127421d9b32SPeter Brune #undef __FUNCT__
128421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
129421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
130421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
131d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
132421d9b32SPeter Brune 
133421d9b32SPeter Brune   PetscFunctionBegin;
134421d9b32SPeter Brune   if (level > top_level)
135421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
136421d9b32SPeter Brune   /* get to the correct level */
137d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
138421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
139421d9b32SPeter Brune   }
140421d9b32SPeter Brune   if (fas->level != level)
141421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
142421d9b32SPeter Brune   fas->restrct = mat;
143421d9b32SPeter Brune   PetscFunctionReturn(0);
144421d9b32SPeter Brune }
145421d9b32SPeter Brune 
146421d9b32SPeter Brune #undef __FUNCT__
147421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
148421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
149421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
150d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
151421d9b32SPeter Brune 
152421d9b32SPeter Brune   PetscFunctionBegin;
153421d9b32SPeter Brune   if (level > top_level)
154421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
155421d9b32SPeter Brune   /* get to the correct level */
156d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
157421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
158421d9b32SPeter Brune   }
159421d9b32SPeter Brune   if (fas->level != level)
160421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
161421d9b32SPeter Brune   fas->rscale = rscale;
162421d9b32SPeter Brune   PetscFunctionReturn(0);
163421d9b32SPeter Brune }
164421d9b32SPeter Brune 
165421d9b32SPeter Brune #undef __FUNCT__
166421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
167421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
168421d9b32SPeter Brune {
16977df8cc4SPeter Brune   PetscErrorCode ierr = 0;
170421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
171421d9b32SPeter Brune 
172421d9b32SPeter Brune   PetscFunctionBegin;
173421d9b32SPeter Brune   /* destroy local data created in SNESSetup_FAS */
174421d9b32SPeter Brune   if (fas->presmooth)    ierr = SNESDestroy(&fas->presmooth);CHKERRQ(ierr);
175421d9b32SPeter Brune   if (fas->postsmooth)   ierr = SNESDestroy(&fas->postsmooth);CHKERRQ(ierr);
176421d9b32SPeter Brune   if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
177421d9b32SPeter Brune   if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
178421d9b32SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
179421d9b32SPeter Brune 
180*fa9694d7SPeter Brune 
181421d9b32SPeter Brune   /* recurse -- reset should NOT destroy the structures -- destroy should destroy the structures recursively */
182421d9b32SPeter Brune   if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr);
183*fa9694d7SPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
184421d9b32SPeter Brune   PetscFunctionReturn(0);
185421d9b32SPeter Brune }
186421d9b32SPeter Brune 
187421d9b32SPeter Brune #undef __FUNCT__
188421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
189421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
190421d9b32SPeter Brune {
191421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
192421d9b32SPeter Brune   PetscErrorCode ierr;
193421d9b32SPeter Brune 
194421d9b32SPeter Brune   PetscFunctionBegin;
195421d9b32SPeter Brune   /* recursively resets and then destroys */
196421d9b32SPeter Brune   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
197421d9b32SPeter Brune   if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
198421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
199421d9b32SPeter Brune   PetscFunctionReturn(0);
200421d9b32SPeter Brune }
201421d9b32SPeter Brune 
202421d9b32SPeter Brune #undef __FUNCT__
203421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
204421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
205421d9b32SPeter Brune {
206421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
207421d9b32SPeter Brune   PetscErrorCode ierr;
208421d9b32SPeter Brune 
209421d9b32SPeter Brune   PetscFunctionBegin;
210*fa9694d7SPeter Brune   ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* the four work vectors are used to transfer stuff BACK */
211421d9b32SPeter Brune   /* gets the solver ready for solution */
212421d9b32SPeter Brune   if (snes->dm) {
213421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
214421d9b32SPeter Brune     if (fas->next) {
215421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
216421d9b32SPeter Brune       if (!fas->interpolate) {
217421d9b32SPeter Brune         if (!fas->next->dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "All levels must be assigned a DM");
218421d9b32SPeter Brune         /* set the interpolation and restriction from the DM */
219421d9b32SPeter Brune         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
220421d9b32SPeter Brune         fas->restrct = fas->interpolate;
221421d9b32SPeter Brune       }
222421d9b32SPeter Brune       /* TODO LATER: Preconditioner setup goes here */
223421d9b32SPeter Brune     }
224421d9b32SPeter Brune   }
225*fa9694d7SPeter Brune   if (fas->next) {
226*fa9694d7SPeter Brune     /* gotta set up the solution vector for this to work */
227*fa9694d7SPeter Brune     ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);
228*fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
229*fa9694d7SPeter Brune   }
230*fa9694d7SPeter Brune   /* got to set them all up at once */
231421d9b32SPeter Brune   PetscFunctionReturn(0);
232421d9b32SPeter Brune }
233421d9b32SPeter Brune 
234421d9b32SPeter Brune #undef __FUNCT__
235421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
236421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
237421d9b32SPeter Brune {
238421d9b32SPeter Brune 
239421d9b32SPeter Brune   PetscErrorCode ierr;
240421d9b32SPeter Brune 
241421d9b32SPeter Brune   PetscFunctionBegin;
242421d9b32SPeter Brune   /* types for the pre and postsmoothers */
243421d9b32SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options");CHKERRQ(ierr);
244421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
245421d9b32SPeter Brune   PetscFunctionReturn(0);
246421d9b32SPeter Brune }
247421d9b32SPeter Brune 
248421d9b32SPeter Brune #undef __FUNCT__
249421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
250421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
251421d9b32SPeter Brune {
252421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
253421d9b32SPeter Brune   PetscBool      iascii;
254421d9b32SPeter Brune   PetscErrorCode ierr;
255421d9b32SPeter Brune   PetscInt levels = fas->levels;
256421d9b32SPeter Brune   PetscInt i;
257421d9b32SPeter Brune 
258421d9b32SPeter Brune   PetscFunctionBegin;
259421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
260421d9b32SPeter Brune   if (iascii) {
261421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
262421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
263421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
264421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
265421d9b32SPeter Brune       if (fas->presmooth) {
266421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
267421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
268421d9b32SPeter Brune         ierr = SNESView(fas->presmooth, viewer);CHKERRQ(ierr);
269421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
270421d9b32SPeter Brune       } else {
271421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
272421d9b32SPeter Brune       }
273421d9b32SPeter Brune       if (fas->postsmooth) {
274421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
275421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
276421d9b32SPeter Brune         ierr = SNESView(fas->postsmooth, viewer);CHKERRQ(ierr);
277421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
278421d9b32SPeter Brune       } else {
279421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
280421d9b32SPeter Brune       }
281421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
282421d9b32SPeter Brune     }
283421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
284421d9b32SPeter Brune   } else {
285421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
286421d9b32SPeter Brune   }
287421d9b32SPeter Brune   PetscFunctionReturn(0);
288421d9b32SPeter Brune }
289421d9b32SPeter Brune 
290421d9b32SPeter Brune /*
291*fa9694d7SPeter Brune 
292*fa9694d7SPeter Brune Defines the FAS cycle as:
293*fa9694d7SPeter Brune 
294*fa9694d7SPeter Brune fine problem: F(x) = 0
295*fa9694d7SPeter Brune coarse problem: F^c(x) = b^c
296*fa9694d7SPeter Brune 
297*fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
298*fa9694d7SPeter Brune 
299*fa9694d7SPeter Brune correction:
300*fa9694d7SPeter Brune 
301*fa9694d7SPeter Brune x = x + I(x^c - Rx)
302*fa9694d7SPeter Brune 
303421d9b32SPeter Brune  */
304*fa9694d7SPeter Brune 
305421d9b32SPeter Brune #undef __FUNCT__
306421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private"
307*fa9694d7SPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec F, Vec X) {
308*fa9694d7SPeter Brune 
309*fa9694d7SPeter Brune   PetscErrorCode ierr;
310*fa9694d7SPeter Brune   Vec X_c, Xo_c, F_c, B_c;
311*fa9694d7SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
312421d9b32SPeter Brune 
313421d9b32SPeter Brune   PetscFunctionBegin;
314*fa9694d7SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
315421d9b32SPeter Brune 
316*fa9694d7SPeter Brune   if (fas->next) {
317*fa9694d7SPeter Brune     X_c  = fas->next->work[0];
318*fa9694d7SPeter Brune     Xo_c = fas->next->work[1];
319*fa9694d7SPeter Brune     F_c  = fas->next->work[2];
320*fa9694d7SPeter Brune     B_c  = fas->next->work[3];
321*fa9694d7SPeter Brune     /* project to coarse */
322*fa9694d7SPeter Brune     ierr = MatRestrict(fas->restrct, X, X_c);CHKERRQ(ierr);
323*fa9694d7SPeter Brune     ierr = MatRestrict(fas->restrct, F, F_c);CHKERRQ(ierr);
324*fa9694d7SPeter Brune     ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
325*fa9694d7SPeter Brune     ierr = VecPointwiseMult(F_c,  fas->rscale, F_c);CHKERRQ(ierr);
326*fa9694d7SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = F^c(Rx) - RF(x) */
327*fa9694d7SPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, B_c);CHKERRQ(ierr);   /* B_c = F(X_c) */
328*fa9694d7SPeter Brune     ierr = VecAXPY(B_c, -1.0, F_c);CHKERRQ(ierr);                     /* B_c = F(X_c) - F_C */
329*fa9694d7SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
330*fa9694d7SPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
331*fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
332*fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
333*fa9694d7SPeter Brune     ierr = MatInterpolate(fas->interpolate, X_c, F);CHKERRQ(ierr);
334*fa9694d7SPeter Brune     ierr = VecAXPY(X, -1.0, F);CHKERRQ(ierr);
335*fa9694d7SPeter Brune   }
336*fa9694d7SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
337*fa9694d7SPeter Brune   /* post-smooth -- just update using the post-smoother */
338*fa9694d7SPeter Brune 
339*fa9694d7SPeter Brune   PetscFunctionReturn(0);
340421d9b32SPeter Brune }
341421d9b32SPeter Brune 
342421d9b32SPeter Brune #undef __FUNCT__
343421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
344421d9b32SPeter Brune 
345421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
346421d9b32SPeter Brune {
347*fa9694d7SPeter Brune   PetscErrorCode ierr;
348*fa9694d7SPeter Brune   Vec X, F;
349421d9b32SPeter Brune   PetscFunctionBegin;
350*fa9694d7SPeter Brune   X = snes->vec_sol;
351*fa9694d7SPeter Brune   F = snes->vec_func;
352*fa9694d7SPeter Brune   ierr = FASCycle_Private(snes, snes->vec_func, snes->vec_sol);CHKERRQ(ierr);
353421d9b32SPeter Brune   PetscFunctionReturn(0);
354421d9b32SPeter Brune }
355