xref: /petsc/src/snes/impls/fas/fas.c (revision 794bee337f0d5af71ccfe93f8f7f88922f7fb150)
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;
52421d9b32SPeter Brune 
53421d9b32SPeter Brune   PetscFunctionReturn(0);
54421d9b32SPeter Brune }
55421d9b32SPeter Brune EXTERN_C_END
56421d9b32SPeter Brune 
57421d9b32SPeter Brune #undef __FUNCT__
58421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
59421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
60421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
61421d9b32SPeter Brune   PetscFunctionBegin;
62ee1fd11aSPeter Brune   *levels = fas->levels;
63421d9b32SPeter Brune   PetscFunctionReturn(0);
64421d9b32SPeter Brune }
65421d9b32SPeter Brune 
66421d9b32SPeter Brune #undef __FUNCT__
67646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles"
68646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
69646217ecSPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
70646217ecSPeter Brune   PetscErrorCode ierr;
71646217ecSPeter Brune   PetscFunctionBegin;
72646217ecSPeter Brune   fas->n_cycles = cycles;
73646217ecSPeter Brune   if (fas->next) {
74646217ecSPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
75646217ecSPeter Brune   }
76646217ecSPeter Brune   PetscFunctionReturn(0);
77646217ecSPeter Brune }
78646217ecSPeter Brune 
79646217ecSPeter Brune 
80646217ecSPeter Brune #undef __FUNCT__
81646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel"
82646217ecSPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, PetscInt cycles) {
83646217ecSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
84646217ecSPeter Brune   PetscInt top_level = fas->level,i;
85646217ecSPeter Brune 
86646217ecSPeter Brune   PetscFunctionBegin;
87646217ecSPeter Brune   if (level > top_level)
88646217ecSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
89646217ecSPeter Brune   /* get to the correct level */
90646217ecSPeter Brune   for (i = fas->level; i > level; i--) {
91646217ecSPeter Brune     fas = (SNES_FAS *)fas->next->data;
92646217ecSPeter Brune   }
93646217ecSPeter Brune   if (fas->level != level)
94646217ecSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
95646217ecSPeter Brune   fas->n_cycles = cycles;
96646217ecSPeter Brune   PetscFunctionReturn(0);
97646217ecSPeter Brune }
98646217ecSPeter Brune 
99646217ecSPeter Brune 
100646217ecSPeter Brune #undef __FUNCT__
101421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
102421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
103421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
104421d9b32SPeter Brune   PetscInt levels = fas->level;
105421d9b32SPeter Brune   PetscInt i;
106421d9b32SPeter Brune   PetscFunctionBegin;
107421d9b32SPeter Brune   *lsnes = snes;
108421d9b32SPeter Brune   if (fas->level < level) {
109421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
110421d9b32SPeter Brune   }
111421d9b32SPeter Brune   if (level > levels - 1) {
112421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
113421d9b32SPeter Brune   }
114421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
115421d9b32SPeter Brune     *lsnes = fas->next;
116421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
117421d9b32SPeter Brune   }
118421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
119421d9b32SPeter Brune   PetscFunctionReturn(0);
120421d9b32SPeter Brune }
121421d9b32SPeter Brune 
122421d9b32SPeter Brune #undef __FUNCT__
123421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
124421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
125421d9b32SPeter Brune   PetscErrorCode ierr;
126421d9b32SPeter Brune   PetscInt i;
127421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
128421d9b32SPeter Brune   MPI_Comm comm;
129421d9b32SPeter Brune   PetscFunctionBegin;
130ee1fd11aSPeter Brune   comm = ((PetscObject)snes)->comm;
131ee1fd11aSPeter Brune   if (levels == fas->levels) {
132ee1fd11aSPeter Brune     if (!comms)
133ee1fd11aSPeter Brune       PetscFunctionReturn(0);
134ee1fd11aSPeter Brune   }
135421d9b32SPeter Brune   /* user has changed the number of levels; reset */
136421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
137421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
138421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
139ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
140421d9b32SPeter Brune   /* setup the finest level */
141421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
142421d9b32SPeter Brune     if (comms) comm = comms[i];
143421d9b32SPeter Brune     fas->level = i;
144421d9b32SPeter Brune     fas->levels = levels;
145ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
146421d9b32SPeter Brune     if (i > 0) {
147421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
1484a6b3fb8SBarry Smith       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
149421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
150794bee33SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
151421d9b32SPeter Brune       fas = (SNES_FAS *)fas->next->data;
152421d9b32SPeter Brune     }
153421d9b32SPeter Brune   }
154421d9b32SPeter Brune   PetscFunctionReturn(0);
155421d9b32SPeter Brune }
156421d9b32SPeter Brune 
157421d9b32SPeter Brune #undef __FUNCT__
158421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
159421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
160421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
161d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
162421d9b32SPeter Brune 
163421d9b32SPeter Brune   PetscFunctionBegin;
164421d9b32SPeter Brune   if (level > top_level)
165421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
166421d9b32SPeter Brune   /* get to the correct level */
167d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
168421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
169421d9b32SPeter Brune   }
170421d9b32SPeter Brune   if (fas->level != level)
171421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
172421d9b32SPeter Brune   fas->interpolate = mat;
173421d9b32SPeter Brune   PetscFunctionReturn(0);
174421d9b32SPeter Brune }
175421d9b32SPeter Brune 
176421d9b32SPeter Brune #undef __FUNCT__
177421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
178421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
179421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
180d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
181421d9b32SPeter Brune 
182421d9b32SPeter Brune   PetscFunctionBegin;
183421d9b32SPeter Brune   if (level > top_level)
184421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
185421d9b32SPeter Brune   /* get to the correct level */
186d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
187421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
188421d9b32SPeter Brune   }
189421d9b32SPeter Brune   if (fas->level != level)
190421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
191421d9b32SPeter Brune   fas->restrct = mat;
192421d9b32SPeter Brune   PetscFunctionReturn(0);
193421d9b32SPeter Brune }
194421d9b32SPeter Brune 
195efe1f98aSPeter Brune 
196421d9b32SPeter Brune #undef __FUNCT__
197421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
198421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
199421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
200d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
201421d9b32SPeter Brune 
202421d9b32SPeter Brune   PetscFunctionBegin;
203421d9b32SPeter Brune   if (level > top_level)
204421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
205421d9b32SPeter Brune   /* get to the correct level */
206d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
207421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
208421d9b32SPeter Brune   }
209421d9b32SPeter Brune   if (fas->level != level)
210421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
211421d9b32SPeter Brune   fas->rscale = rscale;
212421d9b32SPeter Brune   PetscFunctionReturn(0);
213421d9b32SPeter Brune }
214421d9b32SPeter Brune 
215efe1f98aSPeter Brune 
216efe1f98aSPeter Brune #undef __FUNCT__
217efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection"
218efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
219efe1f98aSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
220efe1f98aSPeter Brune   PetscInt top_level = fas->level,i;
221efe1f98aSPeter Brune 
222efe1f98aSPeter Brune   PetscFunctionBegin;
223efe1f98aSPeter Brune   if (level > top_level)
224efe1f98aSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
225efe1f98aSPeter Brune   /* get to the correct level */
226efe1f98aSPeter Brune   for (i = fas->level; i > level; i--) {
227efe1f98aSPeter Brune     fas = (SNES_FAS *)fas->next->data;
228efe1f98aSPeter Brune   }
229efe1f98aSPeter Brune   if (fas->level != level)
230efe1f98aSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
231efe1f98aSPeter Brune   fas->inject = mat;
232efe1f98aSPeter Brune   PetscFunctionReturn(0);
233efe1f98aSPeter Brune }
234efe1f98aSPeter Brune 
235421d9b32SPeter Brune #undef __FUNCT__
236421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
237421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
238421d9b32SPeter Brune {
23977df8cc4SPeter Brune   PetscErrorCode ierr = 0;
240421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
241421d9b32SPeter Brune 
242421d9b32SPeter Brune   PetscFunctionBegin;
243421d9b32SPeter Brune   /* destroy local data created in SNESSetup_FAS */
24451e86f29SPeter Brune #if 0
245293a7e31SPeter Brune   /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */
24651e86f29SPeter Brune #endif
247ee1fd11aSPeter Brune   if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr);
24851e86f29SPeter Brune #if 0
24951e86f29SPeter Brune #endif
250421d9b32SPeter Brune   PetscFunctionReturn(0);
251421d9b32SPeter Brune }
252421d9b32SPeter Brune 
253421d9b32SPeter Brune #undef __FUNCT__
254421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
255421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
256421d9b32SPeter Brune {
257421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
258421d9b32SPeter Brune   PetscErrorCode ierr;
259421d9b32SPeter Brune 
260421d9b32SPeter Brune   PetscFunctionBegin;
261421d9b32SPeter Brune   /* recursively resets and then destroys */
262421d9b32SPeter Brune   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
26351e86f29SPeter Brune   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
26451e86f29SPeter Brune   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
265efe1f98aSPeter Brune   if (fas->inject) {
266efe1f98aSPeter Brune     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
267efe1f98aSPeter Brune   }
26851e86f29SPeter Brune   if (fas->interpolate == fas->restrct) {
26951e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
27051e86f29SPeter Brune     fas->restrct = PETSC_NULL;
27151e86f29SPeter Brune   } else {
27251e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
27351e86f29SPeter Brune     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
27451e86f29SPeter Brune   }
27551e86f29SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
27651e86f29SPeter Brune   if (snes->work)        ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);
277421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
278421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
279ee1fd11aSPeter Brune 
280421d9b32SPeter Brune   PetscFunctionReturn(0);
281421d9b32SPeter Brune }
282421d9b32SPeter Brune 
283421d9b32SPeter Brune #undef __FUNCT__
284421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
285421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
286421d9b32SPeter Brune {
2871a266240SBarry Smith   SNES_FAS       *fas = (SNES_FAS *) snes->data,*tmp;
288421d9b32SPeter Brune   PetscErrorCode ierr;
289efe1f98aSPeter Brune   VecScatter     injscatter;
290421d9b32SPeter Brune 
291421d9b32SPeter Brune   PetscFunctionBegin;
2921a266240SBarry Smith   /* should call the SNESSetFromOptions() only when approriate */
2931a266240SBarry Smith   tmp = fas;
2941a266240SBarry Smith   while (tmp) {
295*6bed07a3SBarry Smith     if (tmp->upsmooth) {ierr = SNESSetFromOptions(tmp->upsmooth);CHKERRQ(ierr);}
296*6bed07a3SBarry Smith     if (tmp->downsmooth) {ierr = SNESSetFromOptions(tmp->downsmooth);CHKERRQ(ierr);}
2971a266240SBarry Smith     tmp = tmp->next ? (SNES_FAS*) tmp->next->data : 0;
2981a266240SBarry Smith   }
2991a266240SBarry Smith 
3001a266240SBarry Smith   if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */
301421d9b32SPeter Brune   /* gets the solver ready for solution */
302646217ecSPeter Brune   if (fas->next) {
303646217ecSPeter Brune     if (snes->ops->computegs) {
304646217ecSPeter Brune       ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
305646217ecSPeter Brune     }
306646217ecSPeter Brune   }
307421d9b32SPeter Brune   if (snes->dm) {
308421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
309421d9b32SPeter Brune     if (fas->next) {
310421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
311ee78dd50SPeter Brune       if (!fas->next->dm) {
312ee78dd50SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
313ee78dd50SPeter Brune         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
314ee78dd50SPeter Brune       }
315421d9b32SPeter Brune       /* set the interpolation and restriction from the DM */
316ee78dd50SPeter Brune       if (!fas->interpolate) {
317421d9b32SPeter Brune         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
318421d9b32SPeter Brune         fas->restrct = fas->interpolate;
319421d9b32SPeter Brune       }
320efe1f98aSPeter Brune       /* set the injection from the DM */
321efe1f98aSPeter Brune       if (!fas->inject) {
322efe1f98aSPeter Brune         ierr = DMGetInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
323efe1f98aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
324efe1f98aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
325efe1f98aSPeter Brune       }
326421d9b32SPeter Brune     }
327ee78dd50SPeter Brune     /* set the DMs of the pre and post-smoothers here */
328*6bed07a3SBarry Smith     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
329*6bed07a3SBarry Smith     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
330421d9b32SPeter Brune   }
331fa9694d7SPeter Brune   if (fas->next) {
332fa9694d7SPeter Brune     /* gotta set up the solution vector for this to work */
333*6bed07a3SBarry Smith     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
334fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
335fa9694d7SPeter Brune   }
336fa9694d7SPeter Brune   /* got to set them all up at once */
337421d9b32SPeter Brune   PetscFunctionReturn(0);
338421d9b32SPeter Brune }
339421d9b32SPeter Brune 
340421d9b32SPeter Brune #undef __FUNCT__
341421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
342421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
343421d9b32SPeter Brune {
344ee78dd50SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
345ee78dd50SPeter Brune   PetscInt levels = 1;
346ee78dd50SPeter Brune   PetscBool flg, monflg;
347421d9b32SPeter Brune   PetscErrorCode ierr;
348ee78dd50SPeter Brune   const char * def_smooth = SNESNRICHARDSON;
349ee78dd50SPeter Brune   char pre_type[256];
350ee78dd50SPeter Brune   char post_type[256];
351ee78dd50SPeter Brune   char                    monfilename[PETSC_MAX_PATH_LEN];
352421d9b32SPeter Brune 
353421d9b32SPeter Brune   PetscFunctionBegin;
354c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
355ee78dd50SPeter Brune 
356ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
357ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
358ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
359c732cbdbSBarry Smith     if (!flg && snes->dm) {
360c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
361c732cbdbSBarry Smith       levels++;
362c732cbdbSBarry Smith     }
363ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
364ee78dd50SPeter Brune   }
365ee78dd50SPeter Brune 
366ee78dd50SPeter Brune   /* type of pre and/or post smoothers -- set both at once */
367ee78dd50SPeter Brune   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
368ee78dd50SPeter Brune   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
369ee78dd50SPeter Brune   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
370ee78dd50SPeter Brune   if (flg) {
371ee78dd50SPeter Brune     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
372ee78dd50SPeter Brune   } else {
373ee78dd50SPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
374ee78dd50SPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
375ee78dd50SPeter Brune   }
376ee78dd50SPeter Brune 
377ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
378ee78dd50SPeter Brune   ierr = PetscOptionsInt("-snes_fas_up_it","Number of upsmoother iterations","PCMGSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
379ee78dd50SPeter Brune   ierr = PetscOptionsInt("-snes_fas_down_it","Number of downsmoother iterations","PCMGSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
380ee78dd50SPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","PCMGSetNumberSmoothUp",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
381ee78dd50SPeter Brune 
382646217ecSPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
383ee78dd50SPeter Brune 
384ee78dd50SPeter Brune   /* other options for the coarsest level */
385ee78dd50SPeter Brune   if (fas->level == 0) {
386ee78dd50SPeter Brune     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&flg);CHKERRQ(ierr);
387ee78dd50SPeter Brune     if (flg) {
388ee78dd50SPeter Brune       ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
389ee78dd50SPeter Brune     } else {
390ee78dd50SPeter Brune       ierr = PetscOptionsList("-snes_fas_coarse_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&flg);CHKERRQ(ierr);
391ee78dd50SPeter Brune       ierr = PetscOptionsList("-snes_fas_coarse_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&flg);CHKERRQ(ierr);
392ee78dd50SPeter Brune     }
393ee78dd50SPeter Brune   }
394ee78dd50SPeter Brune 
395421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
396ee78dd50SPeter Brune   /* setup from the determined types if the smoothers don't exist */
397ee78dd50SPeter Brune   if (!fas->upsmooth) {
39867339d5cSBarry Smith     const char     *prefix;
39967339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
400ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
40167339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
4021a266240SBarry Smith     if (fas->level || (fas->levels == 1)) {
40367339d5cSBarry Smith       ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_");CHKERRQ(ierr);
40467339d5cSBarry Smith     } else {
40567339d5cSBarry Smith       ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_coarse_");CHKERRQ(ierr);
40667339d5cSBarry Smith     }
407293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
408ee78dd50SPeter Brune     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
40969bed9afSBarry Smith     if (snes->ops->computefunction) {
41069bed9afSBarry Smith       ierr = SNESSetFunction(fas->upsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
41169bed9afSBarry Smith     }
412ee78dd50SPeter Brune   }
4131a266240SBarry Smith   if (fas->upsmooth) {
4141a266240SBarry Smith     ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
4151a266240SBarry Smith   }
4161a266240SBarry Smith 
417293a7e31SPeter Brune   if (!fas->downsmooth && fas->level != 0) {
41867339d5cSBarry Smith     const char     *prefix;
41967339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
420ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
42167339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
42267339d5cSBarry Smith     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_");CHKERRQ(ierr);
423293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
424ee78dd50SPeter Brune     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
42569bed9afSBarry Smith     if (snes->ops->computefunction) {
42669bed9afSBarry Smith       ierr = SNESSetFunction(fas->downsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
42769bed9afSBarry Smith     }
428ee78dd50SPeter Brune   }
4291a266240SBarry Smith   if (fas->downsmooth) {
4301a266240SBarry Smith     ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
4311a266240SBarry Smith   }
432ee78dd50SPeter Brune 
433ee78dd50SPeter Brune   if (monflg) {
434646217ecSPeter Brune     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
435794bee33SPeter Brune 
436794bee33SPeter Brune     /* set the monitors for the upsmoother and downsmoother */
437293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
438293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
439ee78dd50SPeter Brune   }
440ee78dd50SPeter Brune 
441ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
4421a266240SBarry Smith   if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);}
443421d9b32SPeter Brune   PetscFunctionReturn(0);
444421d9b32SPeter Brune }
445421d9b32SPeter Brune 
446421d9b32SPeter Brune #undef __FUNCT__
447421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
448421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
449421d9b32SPeter Brune {
450421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
451421d9b32SPeter Brune   PetscBool      iascii;
452421d9b32SPeter Brune   PetscErrorCode ierr;
453421d9b32SPeter Brune   PetscInt levels = fas->levels;
454421d9b32SPeter Brune   PetscInt i;
455421d9b32SPeter Brune 
456421d9b32SPeter Brune   PetscFunctionBegin;
457421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
458421d9b32SPeter Brune   if (iascii) {
459421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
460421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
461421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
462421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
463ee78dd50SPeter Brune       if (fas->upsmooth) {
464ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
465421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
466ee78dd50SPeter Brune         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
467421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
468421d9b32SPeter Brune       } else {
469ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
470421d9b32SPeter Brune       }
471ee78dd50SPeter Brune       if (fas->downsmooth) {
472ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
473421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
474ee78dd50SPeter Brune         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
475421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
476421d9b32SPeter Brune       } else {
477ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
478421d9b32SPeter Brune       }
479421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
480421d9b32SPeter Brune     }
481421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
482421d9b32SPeter Brune   } else {
483421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
484421d9b32SPeter Brune   }
485421d9b32SPeter Brune   PetscFunctionReturn(0);
486421d9b32SPeter Brune }
487421d9b32SPeter Brune 
488421d9b32SPeter Brune /*
489fa9694d7SPeter Brune 
490fa9694d7SPeter Brune Defines the FAS cycle as:
491fa9694d7SPeter Brune 
492fa9694d7SPeter Brune fine problem: F(x) = 0
493fa9694d7SPeter Brune coarse problem: F^c(x) = b^c
494fa9694d7SPeter Brune 
495fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
496fa9694d7SPeter Brune 
497fa9694d7SPeter Brune correction:
498fa9694d7SPeter Brune 
499fa9694d7SPeter Brune x = x + I(x^c - Rx)
500fa9694d7SPeter Brune 
501421d9b32SPeter Brune  */
502fa9694d7SPeter Brune 
503421d9b32SPeter Brune #undef __FUNCT__
504421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private"
505646217ecSPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec X) {
506fa9694d7SPeter Brune 
507fa9694d7SPeter Brune   PetscErrorCode ierr;
508646217ecSPeter Brune   Vec X_c, Xo_c, F_c, B_c,F,B;
509fa9694d7SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
510646217ecSPeter Brune   PetscInt i, k;
511646217ecSPeter Brune   PetscScalar fnorm;
512421d9b32SPeter Brune 
513421d9b32SPeter Brune   PetscFunctionBegin;
514f5a6d4f9SBarry Smith   F = snes->vec_func;
515646217ecSPeter Brune   B = snes->vec_rhs;
516fa9694d7SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
517646217ecSPeter Brune   if (snes->ops->computegs) {
518794bee33SPeter Brune     if (fas->monitor) {
519794bee33SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
520794bee33SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
521794bee33SPeter Brune       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
522794bee33SPeter Brune       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d Gauss-Seidel Function norm %e\n", 0, fnorm);CHKERRQ(ierr);
523794bee33SPeter Brune       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
524794bee33SPeter Brune     }
525646217ecSPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
526646217ecSPeter Brune       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
527794bee33SPeter Brune       if (snes->monitor) {
528794bee33SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
529794bee33SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
530794bee33SPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
531794bee33SPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d Gauss-Seidel Function norm %e\n", k+1, fnorm);CHKERRQ(ierr);
532794bee33SPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
533794bee33SPeter Brune       }
534646217ecSPeter Brune     }
535646217ecSPeter Brune   } else if (fas->upsmooth) {
536ee78dd50SPeter Brune     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
537c90fad12SPeter Brune   } else if (snes->pc) {
538c90fad12SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
539fe6f9142SPeter Brune   }
540fa9694d7SPeter Brune   if (fas->next) {
541293a7e31SPeter Brune     for (i = 0; i < fas->n_cycles; i++) {
542ee78dd50SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
543794bee33SPeter Brune 
544c90fad12SPeter Brune       X_c  = fas->next->vec_sol;
545293a7e31SPeter Brune       Xo_c = fas->next->work[0];
546c90fad12SPeter Brune       F_c  = fas->next->vec_func;
547293a7e31SPeter Brune       B_c  = fas->next->work[1];
548efe1f98aSPeter Brune 
549efe1f98aSPeter Brune       /* inject the solution */
550efe1f98aSPeter Brune       if (fas->inject) {
551a3cb90a9SPeter Brune         ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
552efe1f98aSPeter Brune       } else {
553a3cb90a9SPeter Brune         ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
554a3cb90a9SPeter Brune         ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
555efe1f98aSPeter Brune       }
556293a7e31SPeter Brune       ierr = VecScale(F, -1.0);CHKERRQ(ierr);
557293a7e31SPeter Brune 
558293a7e31SPeter Brune       /* restrict the defect */
559293a7e31SPeter Brune       ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
560293a7e31SPeter Brune 
561c90fad12SPeter Brune       /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
562ee78dd50SPeter Brune       fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
563c90fad12SPeter Brune       ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
564293a7e31SPeter Brune       ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
565c90fad12SPeter Brune 
566ee78dd50SPeter Brune       /* set initial guess of the coarse problem to the projected fine solution */
567ee78dd50SPeter Brune       ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
568c90fad12SPeter Brune 
569c90fad12SPeter Brune       /* recurse to the next level */
570f5a6d4f9SBarry Smith       fas->next->vec_rhs = B_c;
571646217ecSPeter Brune       ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr);
572f5a6d4f9SBarry Smith       fas->next->vec_rhs = PETSC_NULL;
573ee78dd50SPeter Brune 
574fa9694d7SPeter Brune       /* correct as x <- x + I(x^c - Rx)*/
575fa9694d7SPeter Brune       ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
576ee78dd50SPeter Brune       ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr);
577293a7e31SPeter Brune     }
578fa9694d7SPeter Brune   }
579ee78dd50SPeter Brune     /* down-smooth -- just update using the down-smoother */
580c90fad12SPeter Brune   if (fas->level != 0) {
581646217ecSPeter Brune     if (snes->ops->computegs) {
582794bee33SPeter Brune       if (fas->monitor) {
583794bee33SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
584794bee33SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
585794bee33SPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
586794bee33SPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d Gauss-Seidel Function norm %e\n", 0, fnorm);CHKERRQ(ierr);
587794bee33SPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
588794bee33SPeter Brune       }
589646217ecSPeter Brune       for (k = 0; k < fas->max_down_it; k++) {
590646217ecSPeter Brune         ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
591794bee33SPeter Brune         if (fas->monitor) {
592794bee33SPeter Brune           ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
593794bee33SPeter Brune           ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
594794bee33SPeter Brune           ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
595794bee33SPeter Brune           ierr = PetscViewerASCIIPrintf(fas->monitor, "%d Gauss-Seidel Function norm %e\n", k+1, fnorm);CHKERRQ(ierr);
596794bee33SPeter Brune           ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 1);CHKERRQ(ierr);
597794bee33SPeter Brune         }
598646217ecSPeter Brune       }
599646217ecSPeter Brune     }else if (fas->downsmooth) {
600ee78dd50SPeter Brune       ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
601c90fad12SPeter Brune     } else if (snes->pc) {
602c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
603c90fad12SPeter Brune     }
604fe6f9142SPeter Brune   }
605fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
606fa9694d7SPeter Brune 
607fa9694d7SPeter Brune   PetscFunctionReturn(0);
608421d9b32SPeter Brune }
609421d9b32SPeter Brune 
610421d9b32SPeter Brune #undef __FUNCT__
611421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
612421d9b32SPeter Brune 
613421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
614421d9b32SPeter Brune {
615fa9694d7SPeter Brune   PetscErrorCode ierr;
616fe6f9142SPeter Brune   PetscInt i, maxits;
617f5a6d4f9SBarry Smith   Vec X, B,F;
618fe6f9142SPeter Brune   PetscReal fnorm;
619421d9b32SPeter Brune   PetscFunctionBegin;
620fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
621fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
622fa9694d7SPeter Brune   X = snes->vec_sol;
623fe6f9142SPeter Brune   B = snes->vec_rhs;
624f5a6d4f9SBarry Smith   F = snes->vec_func;
625293a7e31SPeter Brune 
626293a7e31SPeter Brune   /*norm setup */
627fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
628fe6f9142SPeter Brune   snes->iter = 0;
629fe6f9142SPeter Brune   snes->norm = 0.;
630fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
631fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
632fe6f9142SPeter Brune   if (snes->domainerror) {
633fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
634fe6f9142SPeter Brune     PetscFunctionReturn(0);
635fe6f9142SPeter Brune   }
636fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
637fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
638fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
639fe6f9142SPeter Brune   snes->norm = fnorm;
640fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
641fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
642fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
643fe6f9142SPeter Brune 
644fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
645fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
646fe6f9142SPeter Brune   /* test convergence */
647fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
648fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
649fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
650fe6f9142SPeter Brune     /* Call general purpose update function */
651646217ecSPeter Brune 
652fe6f9142SPeter Brune     if (snes->ops->update) {
653fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
654fe6f9142SPeter Brune     }
655646217ecSPeter Brune     ierr = FASCycle_Private(snes, X);CHKERRQ(ierr);
656c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
657c90fad12SPeter Brune     /* Monitor convergence */
658c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
659c90fad12SPeter Brune     snes->iter = i+1;
660c90fad12SPeter Brune     snes->norm = fnorm;
661c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
662c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
663c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
664c90fad12SPeter Brune     /* Test for convergence */
665c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
666c90fad12SPeter Brune     if (snes->reason) break;
667fe6f9142SPeter Brune   }
668fe6f9142SPeter Brune   if (i == maxits) {
669fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
670fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
671fe6f9142SPeter Brune   }
672421d9b32SPeter Brune   PetscFunctionReturn(0);
673421d9b32SPeter Brune }
674