xref: /petsc/src/snes/impls/fas/fas.c (revision c90fad1263365a1c13a4b04adb9a18fc33a0fc39)
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);
102*c90fad12SPeter Brune       ierr = SNESSetOptionsPrefix(fas->next, "fas_");CHKERRQ(ierr);
103*c90fad12SPeter Brune       ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);
104421d9b32SPeter Brune       fas = (SNES_FAS *)fas->next->data;
105421d9b32SPeter Brune     }
106421d9b32SPeter Brune   }
107421d9b32SPeter Brune   PetscFunctionReturn(0);
108421d9b32SPeter Brune }
109421d9b32SPeter Brune 
110421d9b32SPeter Brune #undef __FUNCT__
111421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
112421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
113421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
114d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
115421d9b32SPeter Brune 
116421d9b32SPeter Brune   PetscFunctionBegin;
117421d9b32SPeter Brune   if (level > top_level)
118421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
119421d9b32SPeter Brune   /* get to the correct level */
120d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
121421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
122421d9b32SPeter Brune   }
123421d9b32SPeter Brune   if (fas->level != level)
124421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
125421d9b32SPeter Brune   fas->interpolate = mat;
126421d9b32SPeter Brune   PetscFunctionReturn(0);
127421d9b32SPeter Brune }
128421d9b32SPeter Brune 
129421d9b32SPeter Brune #undef __FUNCT__
130421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
131421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
132421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
133d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
134421d9b32SPeter Brune 
135421d9b32SPeter Brune   PetscFunctionBegin;
136421d9b32SPeter Brune   if (level > top_level)
137421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
138421d9b32SPeter Brune   /* get to the correct level */
139d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
140421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
141421d9b32SPeter Brune   }
142421d9b32SPeter Brune   if (fas->level != level)
143421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
144421d9b32SPeter Brune   fas->restrct = mat;
145421d9b32SPeter Brune   PetscFunctionReturn(0);
146421d9b32SPeter Brune }
147421d9b32SPeter Brune 
148421d9b32SPeter Brune #undef __FUNCT__
149421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
150421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
151421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
152d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
153421d9b32SPeter Brune 
154421d9b32SPeter Brune   PetscFunctionBegin;
155421d9b32SPeter Brune   if (level > top_level)
156421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
157421d9b32SPeter Brune   /* get to the correct level */
158d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
159421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
160421d9b32SPeter Brune   }
161421d9b32SPeter Brune   if (fas->level != level)
162421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
163421d9b32SPeter Brune   fas->rscale = rscale;
164421d9b32SPeter Brune   PetscFunctionReturn(0);
165421d9b32SPeter Brune }
166421d9b32SPeter Brune 
167421d9b32SPeter Brune #undef __FUNCT__
168421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
169421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
170421d9b32SPeter Brune {
17177df8cc4SPeter Brune   PetscErrorCode ierr = 0;
172421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
173421d9b32SPeter Brune 
174421d9b32SPeter Brune   PetscFunctionBegin;
175421d9b32SPeter Brune   /* destroy local data created in SNESSetup_FAS */
176421d9b32SPeter Brune   if (fas->presmooth)    ierr = SNESDestroy(&fas->presmooth);CHKERRQ(ierr);
177421d9b32SPeter Brune   if (fas->postsmooth)   ierr = SNESDestroy(&fas->postsmooth);CHKERRQ(ierr);
178421d9b32SPeter Brune   if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
179421d9b32SPeter Brune   if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
180421d9b32SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
181421d9b32SPeter Brune 
182421d9b32SPeter Brune   /* recurse -- reset should NOT destroy the structures -- destroy should destroy the structures recursively */
183421d9b32SPeter Brune   if (fas->next) ierr = SNESReset(fas->next);CHKERRQ(ierr);
184fa9694d7SPeter Brune   if (snes->work) {ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);}
185421d9b32SPeter Brune   PetscFunctionReturn(0);
186421d9b32SPeter Brune }
187421d9b32SPeter Brune 
188421d9b32SPeter Brune #undef __FUNCT__
189421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
190421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
191421d9b32SPeter Brune {
192421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
193421d9b32SPeter Brune   PetscErrorCode ierr;
194421d9b32SPeter Brune 
195421d9b32SPeter Brune   PetscFunctionBegin;
196421d9b32SPeter Brune   /* recursively resets and then destroys */
197421d9b32SPeter Brune   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
198421d9b32SPeter Brune   if (fas->next) ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
199421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
200421d9b32SPeter Brune   PetscFunctionReturn(0);
201421d9b32SPeter Brune }
202421d9b32SPeter Brune 
203421d9b32SPeter Brune #undef __FUNCT__
204421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
205421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
206421d9b32SPeter Brune {
207421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
208421d9b32SPeter Brune   PetscErrorCode ierr;
209421d9b32SPeter Brune 
210421d9b32SPeter Brune   PetscFunctionBegin;
211fa9694d7SPeter Brune   ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* the four work vectors are used to transfer stuff BACK */
212421d9b32SPeter Brune   /* gets the solver ready for solution */
213421d9b32SPeter Brune   if (snes->dm) {
214421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
215421d9b32SPeter Brune     if (fas->next) {
216421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
217421d9b32SPeter Brune       if (!fas->interpolate) {
218421d9b32SPeter Brune         if (!fas->next->dm) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONGSTATE, "All levels must be assigned a DM");
219421d9b32SPeter Brune         /* set the interpolation and restriction from the DM */
220421d9b32SPeter Brune         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
221421d9b32SPeter Brune         fas->restrct = fas->interpolate;
222421d9b32SPeter Brune       }
223421d9b32SPeter Brune       /* TODO LATER: Preconditioner setup goes here */
224421d9b32SPeter Brune     }
225421d9b32SPeter Brune   }
226fa9694d7SPeter Brune   if (fas->next) {
227fa9694d7SPeter Brune     /* gotta set up the solution vector for this to work */
228fa9694d7SPeter Brune     ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);
229fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
230fa9694d7SPeter Brune   }
231fa9694d7SPeter Brune   /* got to set them all up at once */
232421d9b32SPeter Brune   PetscFunctionReturn(0);
233421d9b32SPeter Brune }
234421d9b32SPeter Brune 
235421d9b32SPeter Brune #undef __FUNCT__
236421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
237421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
238421d9b32SPeter Brune {
239421d9b32SPeter Brune 
240421d9b32SPeter Brune   PetscErrorCode ierr;
241421d9b32SPeter Brune 
242421d9b32SPeter Brune   PetscFunctionBegin;
243421d9b32SPeter Brune   /* types for the pre and postsmoothers */
244*c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
245421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
246421d9b32SPeter Brune   PetscFunctionReturn(0);
247421d9b32SPeter Brune }
248421d9b32SPeter Brune 
249421d9b32SPeter Brune #undef __FUNCT__
250421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
251421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
252421d9b32SPeter Brune {
253421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
254421d9b32SPeter Brune   PetscBool      iascii;
255421d9b32SPeter Brune   PetscErrorCode ierr;
256421d9b32SPeter Brune   PetscInt levels = fas->levels;
257421d9b32SPeter Brune   PetscInt i;
258421d9b32SPeter Brune 
259421d9b32SPeter Brune   PetscFunctionBegin;
260421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
261421d9b32SPeter Brune   if (iascii) {
262421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
263421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
264421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
265421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
266421d9b32SPeter Brune       if (fas->presmooth) {
267421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
268421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
269421d9b32SPeter Brune         ierr = SNESView(fas->presmooth, viewer);CHKERRQ(ierr);
270421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
271421d9b32SPeter Brune       } else {
272421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no pre-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
273421d9b32SPeter Brune       }
274421d9b32SPeter Brune       if (fas->postsmooth) {
275421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
276421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
277421d9b32SPeter Brune         ierr = SNESView(fas->postsmooth, viewer);CHKERRQ(ierr);
278421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
279421d9b32SPeter Brune       } else {
280421d9b32SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no post-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
281421d9b32SPeter Brune       }
282421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
283421d9b32SPeter Brune     }
284421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
285421d9b32SPeter Brune   } else {
286421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
287421d9b32SPeter Brune   }
288421d9b32SPeter Brune   PetscFunctionReturn(0);
289421d9b32SPeter Brune }
290421d9b32SPeter Brune 
291421d9b32SPeter Brune /*
292fa9694d7SPeter Brune 
293fa9694d7SPeter Brune Defines the FAS cycle as:
294fa9694d7SPeter Brune 
295fa9694d7SPeter Brune fine problem: F(x) = 0
296fa9694d7SPeter Brune coarse problem: F^c(x) = b^c
297fa9694d7SPeter Brune 
298fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
299fa9694d7SPeter Brune 
300fa9694d7SPeter Brune correction:
301fa9694d7SPeter Brune 
302fa9694d7SPeter Brune x = x + I(x^c - Rx)
303fa9694d7SPeter Brune 
304421d9b32SPeter Brune  */
305fa9694d7SPeter Brune 
306421d9b32SPeter Brune #undef __FUNCT__
307421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private"
308*c90fad12SPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec B, Vec X, Vec F) {
309fa9694d7SPeter Brune 
310fa9694d7SPeter Brune   PetscErrorCode ierr;
311fa9694d7SPeter Brune   Vec X_c, Xo_c, F_c, B_c;
312fa9694d7SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
313421d9b32SPeter Brune 
314421d9b32SPeter Brune   PetscFunctionBegin;
315fa9694d7SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
316fe6f9142SPeter Brune   if (fas->presmooth) {
317fe6f9142SPeter Brune     ierr = SNESSolve(fas->presmooth, B, X);CHKERRQ(ierr);
318*c90fad12SPeter Brune   } else if (snes->pc) {
319*c90fad12SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
320fe6f9142SPeter Brune   }
321fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
322fa9694d7SPeter Brune   if (fas->next) {
323*c90fad12SPeter Brune     X_c  = fas->next->vec_sol;
324fa9694d7SPeter Brune     Xo_c = fas->next->work[1];
325*c90fad12SPeter Brune     F_c  = fas->next->vec_func;
326*c90fad12SPeter Brune     B_c  = fas->next->work[2];
327fa9694d7SPeter Brune     /* project to coarse */
328*c90fad12SPeter Brune     ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
329fa9694d7SPeter Brune     ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
330*c90fad12SPeter Brune     ierr = MatRestrict(fas->restrct, F, F_c);CHKERRQ(ierr);
331fa9694d7SPeter Brune     ierr = VecPointwiseMult(F_c,  fas->rscale, F_c);CHKERRQ(ierr);
332*c90fad12SPeter Brune     if (B) {
333*c90fad12SPeter Brune       ierr = MatRestrict(fas->restrct, B, B_c);CHKERRQ(ierr);
334*c90fad12SPeter Brune       ierr = VecPointwiseMult(B_c,  fas->rscale, B_c);CHKERRQ(ierr);
335*c90fad12SPeter Brune     } else {
336*c90fad12SPeter Brune       ierr = VecSet(B_c, 0.0);CHKERRQ(ierr);
337*c90fad12SPeter Brune     }
338*c90fad12SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
339*c90fad12SPeter Brune     fas->next->vec_rhs = PETSC_NULL; /*unset the RHS for the next problem so we may evaluate the function rather than the residual */
340*c90fad12SPeter Brune     ierr = VecAXPY(B_c, -1.0, F_c);CHKERRQ(ierr);                     /* B_c = RB + F(X_c) - F_C */
341*c90fad12SPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
342*c90fad12SPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);
343fa9694d7SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
344*c90fad12SPeter Brune 
345*c90fad12SPeter Brune     /* test -- initial residuals with and without the RHS set */
346*c90fad12SPeter Brune     fas->next->vec_rhs = B_c;
347*c90fad12SPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
348*c90fad12SPeter Brune     fas->next->vec_rhs = PETSC_NULL;
349*c90fad12SPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
350*c90fad12SPeter Brune 
351*c90fad12SPeter Brune     /* recurse to the next level */
352*c90fad12SPeter Brune     ierr = FASCycle_Private(fas->next, B_c, X_c, F_c);CHKERRQ(ierr);
353fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
354fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
355fa9694d7SPeter Brune     ierr = MatInterpolate(fas->interpolate, X_c, F);CHKERRQ(ierr);
356*c90fad12SPeter Brune     ierr = VecAXPY(X, 1.0, F);CHKERRQ(ierr);
357fa9694d7SPeter Brune   }
358fa9694d7SPeter Brune   /* post-smooth -- just update using the post-smoother */
359*c90fad12SPeter Brune   if (fas->level != 0) {
360fe6f9142SPeter Brune     if (fas->postsmooth) {
361*c90fad12SPeter Brune       ierr = SNESSolve(fas->postsmooth, B, X);CHKERRQ(ierr);
362*c90fad12SPeter Brune     } else if (snes->pc) {
363*c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
364*c90fad12SPeter Brune     }
365fe6f9142SPeter Brune   }
366fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
367fa9694d7SPeter Brune 
368fa9694d7SPeter Brune   PetscFunctionReturn(0);
369421d9b32SPeter Brune }
370421d9b32SPeter Brune 
371421d9b32SPeter Brune #undef __FUNCT__
372*c90fad12SPeter Brune #define __FUNCT__ "FASInitialGuess_Private"
373*c90fad12SPeter Brune 
374*c90fad12SPeter Brune 
375*c90fad12SPeter Brune #undef __FUNCT__
376421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
377421d9b32SPeter Brune 
378421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
379421d9b32SPeter Brune {
380fa9694d7SPeter Brune   PetscErrorCode ierr;
381fe6f9142SPeter Brune   PetscInt i, maxits;
382fe6f9142SPeter Brune   Vec X, F, B;
383fe6f9142SPeter Brune   PetscReal fnorm;
384421d9b32SPeter Brune   PetscFunctionBegin;
385fe6f9142SPeter Brune   maxits = snes->max_its;	     /* maximum number of iterations */
386fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
387fa9694d7SPeter Brune   X = snes->vec_sol;
388fa9694d7SPeter Brune   F = snes->vec_func;
389fe6f9142SPeter Brune   B = snes->vec_rhs;
390fe6f9142SPeter Brune   /* initial iteration */
391fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
392fe6f9142SPeter Brune   snes->iter = 0;
393fe6f9142SPeter Brune   snes->norm = 0.;
394fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
395fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
396fe6f9142SPeter Brune   if (snes->domainerror) {
397fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
398fe6f9142SPeter Brune     PetscFunctionReturn(0);
399fe6f9142SPeter Brune   }
400fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
401fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
402fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
403fe6f9142SPeter Brune   snes->norm = fnorm;
404fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
405fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
406fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
407fe6f9142SPeter Brune 
408fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
409fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
410fe6f9142SPeter Brune   /* test convergence */
411fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
412fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
413fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
414fe6f9142SPeter Brune     /* Call general purpose update function */
415fe6f9142SPeter Brune     if (snes->ops->update) {
416fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
417fe6f9142SPeter Brune     }
418*c90fad12SPeter Brune     ierr = FASCycle_Private(snes, B, X, F);CHKERRQ(ierr);
419*c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
420*c90fad12SPeter Brune     /* Monitor convergence */
421*c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
422*c90fad12SPeter Brune     snes->iter = i+1;
423*c90fad12SPeter Brune     snes->norm = fnorm;
424*c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
425*c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
426*c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
427*c90fad12SPeter Brune     /* Test for convergence */
428*c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
429*c90fad12SPeter Brune     if (snes->reason) break;
430fe6f9142SPeter Brune   }
431fe6f9142SPeter Brune   if (i == maxits) {
432fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
433fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
434fe6f9142SPeter Brune   }
435421d9b32SPeter Brune   PetscFunctionReturn(0);
436421d9b32SPeter Brune }
437