xref: /petsc/src/snes/impls/fas/fas.c (revision cc05f8837a687b6976921eccf6154b226b9d3488)
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 
36ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
37ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
38ed020824SBarry Smith 
39421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
40421d9b32SPeter Brune   snes->data                = (void*) fas;
41421d9b32SPeter Brune   fas->level                = 0;
42293a7e31SPeter Brune   fas->levels               = 1;
43ee78dd50SPeter Brune   fas->n_cycles             = 1;
44ee78dd50SPeter Brune   fas->max_up_it            = 1;
45ee78dd50SPeter Brune   fas->max_down_it          = 1;
46ee78dd50SPeter Brune   fas->upsmooth             = PETSC_NULL;
47ee78dd50SPeter Brune   fas->downsmooth           = PETSC_NULL;
48421d9b32SPeter Brune   fas->next                 = PETSC_NULL;
49421d9b32SPeter Brune   fas->interpolate          = PETSC_NULL;
50421d9b32SPeter Brune   fas->restrct              = PETSC_NULL;
51efe1f98aSPeter Brune   fas->inject               = PETSC_NULL;
52cc05f883SPeter Brune   fas->monitor              = PETSC_NULL;
53cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
54421d9b32SPeter Brune 
55421d9b32SPeter Brune   PetscFunctionReturn(0);
56421d9b32SPeter Brune }
57421d9b32SPeter Brune EXTERN_C_END
58421d9b32SPeter Brune 
59421d9b32SPeter Brune #undef __FUNCT__
60421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
61421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
62421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
63421d9b32SPeter Brune   PetscFunctionBegin;
64ee1fd11aSPeter Brune   *levels = fas->levels;
65421d9b32SPeter Brune   PetscFunctionReturn(0);
66421d9b32SPeter Brune }
67421d9b32SPeter Brune 
68421d9b32SPeter Brune #undef __FUNCT__
69646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles"
70646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
71646217ecSPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
72646217ecSPeter Brune   PetscErrorCode ierr;
73646217ecSPeter Brune   PetscFunctionBegin;
74646217ecSPeter Brune   fas->n_cycles = cycles;
75646217ecSPeter Brune   if (fas->next) {
76646217ecSPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
77646217ecSPeter Brune   }
78646217ecSPeter Brune   PetscFunctionReturn(0);
79646217ecSPeter Brune }
80646217ecSPeter Brune 
81646217ecSPeter Brune 
82646217ecSPeter Brune #undef __FUNCT__
83646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel"
842d15c84fSPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) {
85646217ecSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
86646217ecSPeter Brune   PetscInt top_level = fas->level,i;
87646217ecSPeter Brune 
88646217ecSPeter Brune   PetscFunctionBegin;
89646217ecSPeter Brune   if (level > top_level)
902d15c84fSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
91646217ecSPeter Brune   /* get to the correct level */
92646217ecSPeter Brune   for (i = fas->level; i > level; i--) {
93646217ecSPeter Brune     fas = (SNES_FAS *)fas->next->data;
94646217ecSPeter Brune   }
95646217ecSPeter Brune   if (fas->level != level)
962d15c84fSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
97646217ecSPeter Brune   fas->n_cycles = cycles;
98646217ecSPeter Brune   PetscFunctionReturn(0);
99646217ecSPeter Brune }
100646217ecSPeter Brune 
101646217ecSPeter Brune 
102646217ecSPeter Brune #undef __FUNCT__
103421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
104421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
105421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
106421d9b32SPeter Brune   PetscInt levels = fas->level;
107421d9b32SPeter Brune   PetscInt i;
108421d9b32SPeter Brune   PetscFunctionBegin;
109421d9b32SPeter Brune   *lsnes = snes;
110421d9b32SPeter Brune   if (fas->level < level) {
111421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
112421d9b32SPeter Brune   }
113421d9b32SPeter Brune   if (level > levels - 1) {
114421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
115421d9b32SPeter Brune   }
116421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
117421d9b32SPeter Brune     *lsnes = fas->next;
118421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
119421d9b32SPeter Brune   }
120421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
121421d9b32SPeter Brune   PetscFunctionReturn(0);
122421d9b32SPeter Brune }
123421d9b32SPeter Brune 
124421d9b32SPeter Brune #undef __FUNCT__
125421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
126421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
127421d9b32SPeter Brune   PetscErrorCode ierr;
128421d9b32SPeter Brune   PetscInt i;
129421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
130421d9b32SPeter Brune   MPI_Comm comm;
131421d9b32SPeter Brune   PetscFunctionBegin;
132ee1fd11aSPeter Brune   comm = ((PetscObject)snes)->comm;
133ee1fd11aSPeter Brune   if (levels == fas->levels) {
134ee1fd11aSPeter Brune     if (!comms)
135ee1fd11aSPeter Brune       PetscFunctionReturn(0);
136ee1fd11aSPeter Brune   }
137421d9b32SPeter Brune   /* user has changed the number of levels; reset */
138421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
139421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
140421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
141ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
142421d9b32SPeter Brune   /* setup the finest level */
143421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
144421d9b32SPeter Brune     if (comms) comm = comms[i];
145421d9b32SPeter Brune     fas->level = i;
146421d9b32SPeter Brune     fas->levels = levels;
147ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
148421d9b32SPeter Brune     if (i > 0) {
149421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
1504a6b3fb8SBarry Smith       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
151421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
152794bee33SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
153cc05f883SPeter Brune       if (snes->ops->computegs) {
154cc05f883SPeter Brune         ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
155cc05f883SPeter Brune       }
156421d9b32SPeter Brune       fas = (SNES_FAS *)fas->next->data;
157421d9b32SPeter Brune     }
158421d9b32SPeter Brune   }
159421d9b32SPeter Brune   PetscFunctionReturn(0);
160421d9b32SPeter Brune }
161421d9b32SPeter Brune 
162421d9b32SPeter Brune #undef __FUNCT__
163421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
164421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
165421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
166d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
167421d9b32SPeter Brune 
168421d9b32SPeter Brune   PetscFunctionBegin;
169421d9b32SPeter Brune   if (level > top_level)
170421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
171421d9b32SPeter Brune   /* get to the correct level */
172d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
173421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
174421d9b32SPeter Brune   }
175421d9b32SPeter Brune   if (fas->level != level)
176421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
177421d9b32SPeter Brune   fas->interpolate = mat;
178421d9b32SPeter Brune   PetscFunctionReturn(0);
179421d9b32SPeter Brune }
180421d9b32SPeter Brune 
181421d9b32SPeter Brune #undef __FUNCT__
182421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
183421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
184421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
185d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
186421d9b32SPeter Brune 
187421d9b32SPeter Brune   PetscFunctionBegin;
188421d9b32SPeter Brune   if (level > top_level)
189421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
190421d9b32SPeter Brune   /* get to the correct level */
191d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
192421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
193421d9b32SPeter Brune   }
194421d9b32SPeter Brune   if (fas->level != level)
195421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
196421d9b32SPeter Brune   fas->restrct = mat;
197421d9b32SPeter Brune   PetscFunctionReturn(0);
198421d9b32SPeter Brune }
199421d9b32SPeter Brune 
200efe1f98aSPeter Brune 
201421d9b32SPeter Brune #undef __FUNCT__
202421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
203421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
204421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
205d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
206421d9b32SPeter Brune 
207421d9b32SPeter Brune   PetscFunctionBegin;
208421d9b32SPeter Brune   if (level > top_level)
209421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
210421d9b32SPeter Brune   /* get to the correct level */
211d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
212421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
213421d9b32SPeter Brune   }
214421d9b32SPeter Brune   if (fas->level != level)
215421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
216421d9b32SPeter Brune   fas->rscale = rscale;
217421d9b32SPeter Brune   PetscFunctionReturn(0);
218421d9b32SPeter Brune }
219421d9b32SPeter Brune 
220efe1f98aSPeter Brune 
221efe1f98aSPeter Brune #undef __FUNCT__
222efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection"
223efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
224efe1f98aSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
225efe1f98aSPeter Brune   PetscInt top_level = fas->level,i;
226efe1f98aSPeter Brune 
227efe1f98aSPeter Brune   PetscFunctionBegin;
228efe1f98aSPeter Brune   if (level > top_level)
229efe1f98aSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
230efe1f98aSPeter Brune   /* get to the correct level */
231efe1f98aSPeter Brune   for (i = fas->level; i > level; i--) {
232efe1f98aSPeter Brune     fas = (SNES_FAS *)fas->next->data;
233efe1f98aSPeter Brune   }
234efe1f98aSPeter Brune   if (fas->level != level)
235efe1f98aSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
236efe1f98aSPeter Brune   fas->inject = mat;
237efe1f98aSPeter Brune   PetscFunctionReturn(0);
238efe1f98aSPeter Brune }
239efe1f98aSPeter Brune 
240421d9b32SPeter Brune #undef __FUNCT__
241421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
242421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
243421d9b32SPeter Brune {
24477df8cc4SPeter Brune   PetscErrorCode ierr = 0;
245421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
246421d9b32SPeter Brune 
247421d9b32SPeter Brune   PetscFunctionBegin;
248421d9b32SPeter Brune   /* destroy local data created in SNESSetup_FAS */
24951e86f29SPeter Brune #if 0
250293a7e31SPeter Brune   /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */
25151e86f29SPeter Brune #endif
252ee1fd11aSPeter Brune   if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr);
25351e86f29SPeter Brune #if 0
25451e86f29SPeter Brune #endif
255421d9b32SPeter Brune   PetscFunctionReturn(0);
256421d9b32SPeter Brune }
257421d9b32SPeter Brune 
258421d9b32SPeter Brune #undef __FUNCT__
259421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
260421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
261421d9b32SPeter Brune {
262421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
263421d9b32SPeter Brune   PetscErrorCode ierr;
264421d9b32SPeter Brune 
265421d9b32SPeter Brune   PetscFunctionBegin;
266421d9b32SPeter Brune   /* recursively resets and then destroys */
267421d9b32SPeter Brune   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
26851e86f29SPeter Brune   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
26951e86f29SPeter Brune   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
270efe1f98aSPeter Brune   if (fas->inject) {
271efe1f98aSPeter Brune     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
272efe1f98aSPeter Brune   }
27351e86f29SPeter Brune   if (fas->interpolate == fas->restrct) {
27451e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
27551e86f29SPeter Brune     fas->restrct = PETSC_NULL;
27651e86f29SPeter Brune   } else {
27751e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
27851e86f29SPeter Brune     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
27951e86f29SPeter Brune   }
28051e86f29SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
28151e86f29SPeter Brune   if (snes->work)        ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);
282421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
283421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
284ee1fd11aSPeter Brune 
285421d9b32SPeter Brune   PetscFunctionReturn(0);
286421d9b32SPeter Brune }
287421d9b32SPeter Brune 
288421d9b32SPeter Brune #undef __FUNCT__
289421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
290421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
291421d9b32SPeter Brune {
2921a266240SBarry Smith   SNES_FAS       *fas = (SNES_FAS *) snes->data,*tmp;
293421d9b32SPeter Brune   PetscErrorCode ierr;
294efe1f98aSPeter Brune   VecScatter     injscatter;
295d1adcc6fSPeter Brune   PetscInt       dm_levels;
296d1adcc6fSPeter Brune 
297421d9b32SPeter Brune 
298421d9b32SPeter Brune   PetscFunctionBegin;
299d1adcc6fSPeter Brune   /* reset to use the DM if appropriate */
300cc05f883SPeter Brune   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
301d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
302d1adcc6fSPeter Brune     dm_levels++;
303cc05f883SPeter Brune     if (dm_levels > fas->levels) {
304d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
305cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
306d1adcc6fSPeter Brune     }
307d1adcc6fSPeter Brune   }
308d1adcc6fSPeter Brune 
309cc05f883SPeter Brune   if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */
310cc05f883SPeter Brune 
311d1adcc6fSPeter Brune   /* should call the SNESSetFromOptions() only when appropriate */
3121a266240SBarry Smith   tmp = fas;
3131a266240SBarry Smith   while (tmp) {
314*6bed07a3SBarry Smith     if (tmp->upsmooth) {ierr = SNESSetFromOptions(tmp->upsmooth);CHKERRQ(ierr);}
315*6bed07a3SBarry Smith     if (tmp->downsmooth) {ierr = SNESSetFromOptions(tmp->downsmooth);CHKERRQ(ierr);}
3161a266240SBarry Smith     tmp = tmp->next ? (SNES_FAS*) tmp->next->data : 0;
3171a266240SBarry Smith   }
3181a266240SBarry Smith 
319421d9b32SPeter Brune   /* gets the solver ready for solution */
320646217ecSPeter Brune   if (fas->next) {
321646217ecSPeter Brune     if (snes->ops->computegs) {
322646217ecSPeter Brune       ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
323646217ecSPeter Brune     }
324646217ecSPeter Brune   }
325421d9b32SPeter Brune   if (snes->dm) {
326421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
327421d9b32SPeter Brune     if (fas->next) {
328421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
329ee78dd50SPeter Brune       if (!fas->next->dm) {
330ee78dd50SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
331ee78dd50SPeter Brune         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
332ee78dd50SPeter Brune       }
333421d9b32SPeter Brune       /* set the interpolation and restriction from the DM */
334ee78dd50SPeter Brune       if (!fas->interpolate) {
335421d9b32SPeter Brune         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
336421d9b32SPeter Brune         fas->restrct = fas->interpolate;
337421d9b32SPeter Brune       }
338efe1f98aSPeter Brune       /* set the injection from the DM */
339efe1f98aSPeter Brune       if (!fas->inject) {
340efe1f98aSPeter Brune         ierr = DMGetInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
341efe1f98aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
342efe1f98aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
343efe1f98aSPeter Brune       }
344421d9b32SPeter Brune     }
345ee78dd50SPeter Brune     /* set the DMs of the pre and post-smoothers here */
346*6bed07a3SBarry Smith     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
347*6bed07a3SBarry Smith     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
348421d9b32SPeter Brune   }
349cc05f883SPeter Brune 
350fa9694d7SPeter Brune   if (fas->next) {
351fa9694d7SPeter Brune    /* gotta set up the solution vector for this to work */
352*6bed07a3SBarry Smith     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
353fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
354fa9694d7SPeter Brune   }
355fa9694d7SPeter Brune   /* got to set them all up at once */
356421d9b32SPeter Brune   PetscFunctionReturn(0);
357421d9b32SPeter Brune }
358421d9b32SPeter Brune 
359421d9b32SPeter Brune #undef __FUNCT__
360421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
361421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
362421d9b32SPeter Brune {
363ee78dd50SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
364ee78dd50SPeter Brune   PetscInt levels = 1;
365d1adcc6fSPeter Brune   PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg;
366421d9b32SPeter Brune   PetscErrorCode ierr;
367ee78dd50SPeter Brune   const char * def_smooth = SNESNRICHARDSON;
368ee78dd50SPeter Brune   char pre_type[256];
369ee78dd50SPeter Brune   char post_type[256];
370ee78dd50SPeter Brune   char                    monfilename[PETSC_MAX_PATH_LEN];
371421d9b32SPeter Brune 
372421d9b32SPeter Brune   PetscFunctionBegin;
373c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
374ee78dd50SPeter Brune 
375ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
376ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
377ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
378c732cbdbSBarry Smith     if (!flg && snes->dm) {
379c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
380c732cbdbSBarry Smith       levels++;
381d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
382c732cbdbSBarry Smith     }
383ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
384ee78dd50SPeter Brune   }
385ee78dd50SPeter Brune 
386ee78dd50SPeter Brune   /* type of pre and/or post smoothers -- set both at once */
387ee78dd50SPeter Brune   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
388ee78dd50SPeter Brune   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
389d1adcc6fSPeter Brune   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
390d1adcc6fSPeter Brune   if (smoothflg) {
391ee78dd50SPeter Brune     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
392ee78dd50SPeter Brune   } else {
393d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr);
394d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr);
395ee78dd50SPeter Brune   }
396ee78dd50SPeter Brune 
397ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
398d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
399d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
400d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
401ee78dd50SPeter Brune 
402646217ecSPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
403ee78dd50SPeter Brune 
404ee78dd50SPeter Brune   /* other options for the coarsest level */
405ee78dd50SPeter Brune   if (fas->level == 0) {
406d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
407ee78dd50SPeter Brune   }
408ee78dd50SPeter Brune 
409421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
4108cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
4118cc86e31SPeter Brune   if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->ops->computegs)) {
4128cc86e31SPeter Brune     const char     *prefix;
4138cc86e31SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
4148cc86e31SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
4158cc86e31SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
4168cc86e31SPeter Brune     if (fas->level || (fas->levels == 1)) {
4178cc86e31SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_");CHKERRQ(ierr);
4188cc86e31SPeter Brune     } else {
4198cc86e31SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
4208cc86e31SPeter Brune     }
4218cc86e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
4228cc86e31SPeter Brune     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
4238cc86e31SPeter Brune     if (snes->ops->computefunction) {
4248cc86e31SPeter Brune       ierr = SNESSetFunction(fas->downsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
4258cc86e31SPeter Brune     }
4268cc86e31SPeter Brune   }
4278cc86e31SPeter Brune 
4288cc86e31SPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->ops->computegs)) {
42967339d5cSBarry Smith     const char     *prefix;
43067339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
431ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
43267339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
43367339d5cSBarry Smith     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_");CHKERRQ(ierr);
434293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
435ee78dd50SPeter Brune     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
43669bed9afSBarry Smith     if (snes->ops->computefunction) {
43769bed9afSBarry Smith       ierr = SNESSetFunction(fas->upsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
43869bed9afSBarry Smith     }
439ee78dd50SPeter Brune   }
4401a266240SBarry Smith   if (fas->upsmooth) {
4411a266240SBarry Smith     ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
4421a266240SBarry Smith   }
4431a266240SBarry Smith 
4441a266240SBarry Smith   if (fas->downsmooth) {
4451a266240SBarry Smith     ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
4461a266240SBarry Smith   }
447ee78dd50SPeter Brune 
448ee78dd50SPeter Brune   if (monflg) {
449646217ecSPeter Brune     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
450794bee33SPeter Brune     /* set the monitors for the upsmoother and downsmoother */
451293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
452293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
453ee78dd50SPeter Brune   }
454ee78dd50SPeter Brune 
455ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
4561a266240SBarry Smith   if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);}
457421d9b32SPeter Brune   PetscFunctionReturn(0);
458421d9b32SPeter Brune }
459421d9b32SPeter Brune 
460421d9b32SPeter Brune #undef __FUNCT__
461421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
462421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
463421d9b32SPeter Brune {
464421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
465421d9b32SPeter Brune   PetscBool      iascii;
466421d9b32SPeter Brune   PetscErrorCode ierr;
467421d9b32SPeter Brune   PetscInt levels = fas->levels;
468421d9b32SPeter Brune   PetscInt i;
469421d9b32SPeter Brune 
470421d9b32SPeter Brune   PetscFunctionBegin;
471421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
472421d9b32SPeter Brune   if (iascii) {
473421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
474421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
475421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
476421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
477ee78dd50SPeter Brune       if (fas->upsmooth) {
478ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
479421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
480ee78dd50SPeter Brune         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
481421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
482421d9b32SPeter Brune       } else {
483ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
484421d9b32SPeter Brune       }
485ee78dd50SPeter Brune       if (fas->downsmooth) {
486ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
487421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
488ee78dd50SPeter Brune         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
489421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
490421d9b32SPeter Brune       } else {
491ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
492421d9b32SPeter Brune       }
493421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
494421d9b32SPeter Brune     }
495421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
496421d9b32SPeter Brune   } else {
497421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
498421d9b32SPeter Brune   }
499421d9b32SPeter Brune   PetscFunctionReturn(0);
500421d9b32SPeter Brune }
501421d9b32SPeter Brune 
502421d9b32SPeter Brune /*
503fa9694d7SPeter Brune 
504fa9694d7SPeter Brune Defines the FAS cycle as:
505fa9694d7SPeter Brune 
506fa9694d7SPeter Brune fine problem: F(x) = 0
507fa9694d7SPeter Brune coarse problem: F^c(x) = b^c
508fa9694d7SPeter Brune 
509fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
510fa9694d7SPeter Brune 
511fa9694d7SPeter Brune correction:
512fa9694d7SPeter Brune 
513fa9694d7SPeter Brune x = x + I(x^c - Rx)
514fa9694d7SPeter Brune 
515421d9b32SPeter Brune  */
516fa9694d7SPeter Brune 
517421d9b32SPeter Brune #undef __FUNCT__
518421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private"
519646217ecSPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec X) {
520fa9694d7SPeter Brune 
521fa9694d7SPeter Brune   PetscErrorCode ierr;
522646217ecSPeter Brune   Vec X_c, Xo_c, F_c, B_c,F,B;
523fa9694d7SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
524646217ecSPeter Brune   PetscInt i, k;
5252d15c84fSPeter Brune   PetscReal fnorm;
526421d9b32SPeter Brune 
527421d9b32SPeter Brune   PetscFunctionBegin;
528f5a6d4f9SBarry Smith   F = snes->vec_func;
529646217ecSPeter Brune   B = snes->vec_rhs;
530fa9694d7SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
531d1adcc6fSPeter Brune   if (fas->downsmooth) {
532d1adcc6fSPeter Brune     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
533d1adcc6fSPeter Brune   } else if (snes->ops->computegs) {
534794bee33SPeter Brune     if (fas->monitor) {
535794bee33SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
536794bee33SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
537d1adcc6fSPeter Brune       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
538d1adcc6fSPeter Brune       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %e\n", 0, fnorm);CHKERRQ(ierr);
539d1adcc6fSPeter Brune       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
540794bee33SPeter Brune     }
541646217ecSPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
542646217ecSPeter Brune       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
543cc05f883SPeter Brune       if (fas->monitor) {
544794bee33SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
545794bee33SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
546d1adcc6fSPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
547d1adcc6fSPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %e\n", k+1, fnorm);CHKERRQ(ierr);
548d1adcc6fSPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
549794bee33SPeter Brune       }
550646217ecSPeter Brune     }
551c90fad12SPeter Brune   } else if (snes->pc) {
552c90fad12SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
553fe6f9142SPeter Brune   }
554fa9694d7SPeter Brune   if (fas->next) {
555293a7e31SPeter Brune     for (i = 0; i < fas->n_cycles; i++) {
556ee78dd50SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
557794bee33SPeter Brune 
558c90fad12SPeter Brune       X_c  = fas->next->vec_sol;
559293a7e31SPeter Brune       Xo_c = fas->next->work[0];
560c90fad12SPeter Brune       F_c  = fas->next->vec_func;
561293a7e31SPeter Brune       B_c  = fas->next->work[1];
562efe1f98aSPeter Brune 
563efe1f98aSPeter Brune       /* inject the solution */
564efe1f98aSPeter Brune       if (fas->inject) {
565a3cb90a9SPeter Brune         ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
566efe1f98aSPeter Brune       } else {
567a3cb90a9SPeter Brune         ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
568a3cb90a9SPeter Brune         ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
569efe1f98aSPeter Brune       }
570293a7e31SPeter Brune       ierr = VecScale(F, -1.0);CHKERRQ(ierr);
571293a7e31SPeter Brune 
572293a7e31SPeter Brune       /* restrict the defect */
573293a7e31SPeter Brune       ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
574293a7e31SPeter Brune 
575c90fad12SPeter Brune       /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
576ee78dd50SPeter Brune       fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
577c90fad12SPeter Brune       ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
578293a7e31SPeter Brune       ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
579c90fad12SPeter Brune 
580ee78dd50SPeter Brune       /* set initial guess of the coarse problem to the projected fine solution */
581ee78dd50SPeter Brune       ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
582c90fad12SPeter Brune 
583c90fad12SPeter Brune       /* recurse to the next level */
584f5a6d4f9SBarry Smith       fas->next->vec_rhs = B_c;
585646217ecSPeter Brune       ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr);
586f5a6d4f9SBarry Smith       fas->next->vec_rhs = PETSC_NULL;
587ee78dd50SPeter Brune 
588fa9694d7SPeter Brune       /* correct as x <- x + I(x^c - Rx)*/
589fa9694d7SPeter Brune       ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
590ee78dd50SPeter Brune       ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr);
591293a7e31SPeter Brune     }
592fa9694d7SPeter Brune   }
593d1adcc6fSPeter Brune     /* up-smooth -- just update using the down-smoother */
594c90fad12SPeter Brune   if (fas->level != 0) {
595d1adcc6fSPeter Brune     if (fas->upsmooth) {
596d1adcc6fSPeter Brune       ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
597d1adcc6fSPeter Brune     } else if (snes->ops->computegs) {
598794bee33SPeter Brune       if (fas->monitor) {
599794bee33SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
600794bee33SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
601d1adcc6fSPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
602d1adcc6fSPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %e\n", 0, fnorm);CHKERRQ(ierr);
603d1adcc6fSPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
604794bee33SPeter Brune       }
605646217ecSPeter Brune       for (k = 0; k < fas->max_down_it; k++) {
606646217ecSPeter Brune         ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
607794bee33SPeter Brune         if (fas->monitor) {
608794bee33SPeter Brune           ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
609794bee33SPeter Brune           ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
610d1adcc6fSPeter Brune           ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
611d1adcc6fSPeter Brune           ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %e\n", k+1, fnorm);CHKERRQ(ierr);
612d1adcc6fSPeter Brune           ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
613794bee33SPeter Brune         }
614646217ecSPeter Brune       }
615c90fad12SPeter Brune     } else if (snes->pc) {
616c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
617c90fad12SPeter Brune     }
618fe6f9142SPeter Brune   }
619fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
620fa9694d7SPeter Brune 
621fa9694d7SPeter Brune   PetscFunctionReturn(0);
622421d9b32SPeter Brune }
623421d9b32SPeter Brune 
624421d9b32SPeter Brune #undef __FUNCT__
625421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
626421d9b32SPeter Brune 
627421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
628421d9b32SPeter Brune {
629fa9694d7SPeter Brune   PetscErrorCode ierr;
630fe6f9142SPeter Brune   PetscInt i, maxits;
631f5a6d4f9SBarry Smith   Vec X, B,F;
632fe6f9142SPeter Brune   PetscReal fnorm;
633421d9b32SPeter Brune   PetscFunctionBegin;
634fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
635fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
636fa9694d7SPeter Brune   X = snes->vec_sol;
637fe6f9142SPeter Brune   B = snes->vec_rhs;
638f5a6d4f9SBarry Smith   F = snes->vec_func;
639293a7e31SPeter Brune 
640293a7e31SPeter Brune   /*norm setup */
641fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
642fe6f9142SPeter Brune   snes->iter = 0;
643fe6f9142SPeter Brune   snes->norm = 0.;
644fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
645fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
646fe6f9142SPeter Brune   if (snes->domainerror) {
647fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
648fe6f9142SPeter Brune     PetscFunctionReturn(0);
649fe6f9142SPeter Brune   }
650fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
651fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
652fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
653fe6f9142SPeter Brune   snes->norm = fnorm;
654fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
655fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
656fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
657fe6f9142SPeter Brune 
658fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
659fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
660fe6f9142SPeter Brune   /* test convergence */
661fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
662fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
663fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
664fe6f9142SPeter Brune     /* Call general purpose update function */
665646217ecSPeter Brune 
666fe6f9142SPeter Brune     if (snes->ops->update) {
667fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
668fe6f9142SPeter Brune     }
669646217ecSPeter Brune     ierr = FASCycle_Private(snes, X);CHKERRQ(ierr);
670c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
671c90fad12SPeter Brune     /* Monitor convergence */
672c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
673c90fad12SPeter Brune     snes->iter = i+1;
674c90fad12SPeter Brune     snes->norm = fnorm;
675c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
676c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
677c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
678c90fad12SPeter Brune     /* Test for convergence */
679c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
680c90fad12SPeter Brune     if (snes->reason) break;
681fe6f9142SPeter Brune   }
682fe6f9142SPeter Brune   if (i == maxits) {
683fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
684fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
685fe6f9142SPeter Brune   }
686421d9b32SPeter Brune   PetscFunctionReturn(0);
687421d9b32SPeter Brune }
688