xref: /petsc/src/snes/impls/fas/fas.c (revision c79ef2592ddd350a3e361856ed98a4a2544c077c)
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;
54eff52c0eSPeter Brune   fas->useGS                  = PETSC_FALSE;
55421d9b32SPeter Brune   PetscFunctionReturn(0);
56421d9b32SPeter Brune }
57421d9b32SPeter Brune EXTERN_C_END
58421d9b32SPeter Brune 
59421d9b32SPeter Brune #undef __FUNCT__
60421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
61c79ef259SPeter Brune /*@
62c79ef259SPeter Brune    SNESFASGetLevels - Gets the number of levels to use with FAS.
63c79ef259SPeter Brune 
64c79ef259SPeter Brune    Input Parameter:
65c79ef259SPeter Brune .  snes - the preconditioner context
66c79ef259SPeter Brune 
67c79ef259SPeter Brune    Output parameter:
68c79ef259SPeter Brune .  levels - the number of levels
69c79ef259SPeter Brune 
70c79ef259SPeter Brune    Level: advanced
71c79ef259SPeter Brune 
72c79ef259SPeter Brune .keywords: MG, get, levels, multigrid
73c79ef259SPeter Brune 
74c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels()
75c79ef259SPeter Brune @*/
76421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
77421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
78421d9b32SPeter Brune   PetscFunctionBegin;
79ee1fd11aSPeter Brune   *levels = fas->levels;
80421d9b32SPeter Brune   PetscFunctionReturn(0);
81421d9b32SPeter Brune }
82421d9b32SPeter Brune 
83421d9b32SPeter Brune #undef __FUNCT__
84646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles"
85c79ef259SPeter Brune /*@
86c79ef259SPeter Brune    SNESFASSetCycles - Sets the type cycles to use.  Use SNESFASSetCyclesOnLevel() for more
87c79ef259SPeter Brune    complicated cycling.
88c79ef259SPeter Brune 
89c79ef259SPeter Brune    Logically Collective on SNES
90c79ef259SPeter Brune 
91c79ef259SPeter Brune    Input Parameters:
92c79ef259SPeter Brune +  snes   - the multigrid context
93c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
94c79ef259SPeter Brune 
95c79ef259SPeter Brune    Options Database Key:
96c79ef259SPeter Brune $  -snes_fas_cycles 1 or 2
97c79ef259SPeter Brune 
98c79ef259SPeter Brune    Level: advanced
99c79ef259SPeter Brune 
100c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
101c79ef259SPeter Brune 
102c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel()
103c79ef259SPeter Brune @*/
104646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
105646217ecSPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
106646217ecSPeter Brune   PetscErrorCode ierr;
107646217ecSPeter Brune   PetscFunctionBegin;
108646217ecSPeter Brune   fas->n_cycles = cycles;
109646217ecSPeter Brune   if (fas->next) {
110646217ecSPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
111646217ecSPeter Brune   }
112646217ecSPeter Brune   PetscFunctionReturn(0);
113646217ecSPeter Brune }
114646217ecSPeter Brune 
115eff52c0eSPeter Brune #undef __FUNCT__
116c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel"
117c79ef259SPeter Brune /*@
118c79ef259SPeter Brune    SNESFASSetCycles - Sets the type cycles to use.  Use SNESFASSetCyclesOnLevel() for more
119c79ef259SPeter Brune    complicated cycling.
120c79ef259SPeter Brune 
121c79ef259SPeter Brune    Logically Collective on SNES
122c79ef259SPeter Brune 
123c79ef259SPeter Brune    Input Parameters:
124c79ef259SPeter Brune +  snes   - the multigrid context
125c79ef259SPeter Brune .  level  - the level to set the number of cycles on
126c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
127c79ef259SPeter Brune 
128c79ef259SPeter Brune    Level: advanced
129c79ef259SPeter Brune 
130c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
131c79ef259SPeter Brune 
132c79ef259SPeter Brune .seealso: SNESFASSetCycles()
133c79ef259SPeter Brune @*/
134c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) {
135c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
136c79ef259SPeter Brune   PetscInt top_level = fas->level,i;
137c79ef259SPeter Brune 
138c79ef259SPeter Brune   PetscFunctionBegin;
139c79ef259SPeter Brune   if (level > top_level)
140c79ef259SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
141c79ef259SPeter Brune   /* get to the correct level */
142c79ef259SPeter Brune   for (i = fas->level; i > level; i--) {
143c79ef259SPeter Brune     fas = (SNES_FAS *)fas->next->data;
144c79ef259SPeter Brune   }
145c79ef259SPeter Brune   if (fas->level != level)
146c79ef259SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
147c79ef259SPeter Brune   fas->n_cycles = cycles;
148c79ef259SPeter Brune   PetscFunctionReturn(0);
149c79ef259SPeter Brune }
150c79ef259SPeter Brune 
151c79ef259SPeter Brune #undef __FUNCT__
152eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS"
153c79ef259SPeter Brune /*@
154c79ef259SPeter Brune    SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used.
155c79ef259SPeter Brune    Use SNESFASSetGSOnLevel() for more complicated staging of smoothers
156c79ef259SPeter Brune    and nonlinear preconditioners.
157c79ef259SPeter Brune 
158c79ef259SPeter Brune    Logically Collective on SNES
159c79ef259SPeter Brune 
160c79ef259SPeter Brune    Input Parameters:
161c79ef259SPeter Brune +  snes    - the multigrid context
162c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
163c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
164c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
165c79ef259SPeter Brune 
166c79ef259SPeter Brune    Level: advanced
167c79ef259SPeter Brune 
168c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid
169c79ef259SPeter Brune 
170c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel()
171c79ef259SPeter Brune @*/
172eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
173eff52c0eSPeter Brune   PetscErrorCode ierr = 0;
174eff52c0eSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
175eff52c0eSPeter Brune   PetscFunctionBegin;
176eff52c0eSPeter Brune 
177eff52c0eSPeter Brune   /* use or don't use it according to user wishes*/
178eff52c0eSPeter Brune   fas->useGS = use_gs;
179eff52c0eSPeter Brune   if (gsfunc) {
180eff52c0eSPeter Brune     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
181eff52c0eSPeter Brune     /* push the provided GS up the tree */
182eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
183eff52c0eSPeter Brune   } else if (snes->ops->computegs) {
184eff52c0eSPeter Brune     /* assume that the user has set the GS solver at this level */
185eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
186eff52c0eSPeter Brune   } else if (use_gs) {
187eff52c0eSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level);
188eff52c0eSPeter Brune   }
189eff52c0eSPeter Brune   PetscFunctionReturn(0);
190eff52c0eSPeter Brune }
191eff52c0eSPeter Brune 
192eff52c0eSPeter Brune #undef __FUNCT__
193eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel"
194c79ef259SPeter Brune /*@
195c79ef259SPeter Brune    SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level.
196c79ef259SPeter Brune 
197c79ef259SPeter Brune    Logically Collective on SNES
198c79ef259SPeter Brune 
199c79ef259SPeter Brune    Input Parameters:
200c79ef259SPeter Brune +  snes    - the multigrid context
201c79ef259SPeter Brune .  level   - the level to set the nonlinear smoother on
202c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
203c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
204c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
205c79ef259SPeter Brune 
206c79ef259SPeter Brune    Level: advanced
207c79ef259SPeter Brune 
208c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
209c79ef259SPeter Brune 
210c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS()
211c79ef259SPeter Brune @*/
212eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
213eff52c0eSPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
214eff52c0eSPeter Brune   PetscErrorCode ierr;
215eff52c0eSPeter Brune   PetscInt       top_level = fas->level,i;
216eff52c0eSPeter Brune   SNES           cur_snes = snes;
217eff52c0eSPeter Brune   PetscFunctionBegin;
218eff52c0eSPeter Brune   if (level > top_level)
219eff52c0eSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
220eff52c0eSPeter Brune   /* get to the correct level */
221eff52c0eSPeter Brune   for (i = fas->level; i > level; i--) {
222eff52c0eSPeter Brune     fas = (SNES_FAS *)fas->next->data;
223eff52c0eSPeter Brune     cur_snes = fas->next;
224eff52c0eSPeter Brune   }
225eff52c0eSPeter Brune   if (fas->level != level)
226eff52c0eSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
227eff52c0eSPeter Brune   fas->useGS = use_gs;
228eff52c0eSPeter Brune   if (gsfunc) {
229eff52c0eSPeter Brune     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
230eff52c0eSPeter Brune   }
231eff52c0eSPeter Brune   PetscFunctionReturn(0);
232eff52c0eSPeter Brune }
233eff52c0eSPeter Brune 
234646217ecSPeter Brune #undef __FUNCT__
235421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
236c79ef259SPeter Brune /*@
237c79ef259SPeter Brune    SNESFASGetSNES - Gets the SNES corresponding to a particular
238c79ef259SPeter Brune    level of the FAS hierarchy.
239c79ef259SPeter Brune 
240c79ef259SPeter Brune    Input Parameters:
241c79ef259SPeter Brune +  snes    - the multigrid context
242c79ef259SPeter Brune    level   - the level to get
243c79ef259SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
244c79ef259SPeter Brune 
245c79ef259SPeter Brune    Level: advanced
246c79ef259SPeter Brune 
247c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
248c79ef259SPeter Brune 
249c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
250c79ef259SPeter Brune @*/
251421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
252421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
253421d9b32SPeter Brune   PetscInt levels = fas->level;
254421d9b32SPeter Brune   PetscInt i;
255421d9b32SPeter Brune   PetscFunctionBegin;
256421d9b32SPeter Brune   *lsnes = snes;
257421d9b32SPeter Brune   if (fas->level < level) {
258421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
259421d9b32SPeter Brune   }
260421d9b32SPeter Brune   if (level > levels - 1) {
261421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
262421d9b32SPeter Brune   }
263421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
264421d9b32SPeter Brune     *lsnes = fas->next;
265421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
266421d9b32SPeter Brune   }
267421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
268421d9b32SPeter Brune   PetscFunctionReturn(0);
269421d9b32SPeter Brune }
270421d9b32SPeter Brune 
271421d9b32SPeter Brune #undef __FUNCT__
272421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
273c79ef259SPeter Brune /*@
274c79ef259SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
275c79ef259SPeter Brune    Must be called before any other FAS routine.
276c79ef259SPeter Brune 
277c79ef259SPeter Brune    Input Parameters:
278c79ef259SPeter Brune +  snes   - the snes context
279c79ef259SPeter Brune .  levels - the number of levels
280c79ef259SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
281c79ef259SPeter Brune             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
282c79ef259SPeter Brune             Fortran.
283c79ef259SPeter Brune 
284c79ef259SPeter Brune    Level: intermediate
285c79ef259SPeter Brune 
286c79ef259SPeter Brune    Notes:
287c79ef259SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
288c79ef259SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
289c79ef259SPeter Brune 
290c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid
291c79ef259SPeter Brune 
292c79ef259SPeter Brune .seealso: SNESFASGetLevels()
293c79ef259SPeter Brune @*/
294421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
295421d9b32SPeter Brune   PetscErrorCode ierr;
296421d9b32SPeter Brune   PetscInt i;
297421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
298421d9b32SPeter Brune   MPI_Comm comm;
299421d9b32SPeter Brune   PetscFunctionBegin;
300ee1fd11aSPeter Brune    comm = ((PetscObject)snes)->comm;
301ee1fd11aSPeter Brune   if (levels == fas->levels) {
302ee1fd11aSPeter Brune     if (!comms)
303ee1fd11aSPeter Brune       PetscFunctionReturn(0);
304ee1fd11aSPeter Brune   }
305421d9b32SPeter Brune   /* user has changed the number of levels; reset */
306421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
307421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
308421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
309ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
310421d9b32SPeter Brune   /* setup the finest level */
311421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
312421d9b32SPeter Brune     if (comms) comm = comms[i];
313421d9b32SPeter Brune     fas->level = i;
314421d9b32SPeter Brune     fas->levels = levels;
315ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
316421d9b32SPeter Brune     if (i > 0) {
317421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
3184a6b3fb8SBarry Smith       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
319421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
320794bee33SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
321421d9b32SPeter Brune       fas = (SNES_FAS *)fas->next->data;
322421d9b32SPeter Brune     }
323421d9b32SPeter Brune   }
324eff52c0eSPeter Brune   /* by default set the GS smoothers up the chain if one exists on the finest level */
325eff52c0eSPeter Brune   if (snes->ops->computegs)
326eff52c0eSPeter Brune     ierr = SNESFASSetGS(snes, snes->ops->computegs, snes->gsP, fas->useGS);CHKERRQ(ierr);
327421d9b32SPeter Brune   PetscFunctionReturn(0);
328421d9b32SPeter Brune }
329421d9b32SPeter Brune 
330421d9b32SPeter Brune #undef __FUNCT__
331c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp"
332c79ef259SPeter Brune /*@
333c79ef259SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
334c79ef259SPeter Brune    use on all levels.
335c79ef259SPeter Brune 
336c79ef259SPeter Brune    Logically Collective on PC
337c79ef259SPeter Brune 
338c79ef259SPeter Brune    Input Parameters:
339c79ef259SPeter Brune +  snes - the multigrid context
340c79ef259SPeter Brune -  n    - the number of smoothing steps
341c79ef259SPeter Brune 
342c79ef259SPeter Brune    Options Database Key:
343c79ef259SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
344c79ef259SPeter Brune 
345c79ef259SPeter Brune    Level: advanced
346c79ef259SPeter Brune 
347c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
348c79ef259SPeter Brune 
349c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
350c79ef259SPeter Brune @*/
351c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
352c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
353c79ef259SPeter Brune   PetscErrorCode ierr = 0;
354c79ef259SPeter Brune   PetscFunctionBegin;
355c79ef259SPeter Brune   fas->max_up_it = n;
356c79ef259SPeter Brune   if (fas->next) {
357c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
358c79ef259SPeter Brune   }
359c79ef259SPeter Brune   PetscFunctionReturn(0);
360c79ef259SPeter Brune }
361c79ef259SPeter Brune 
362c79ef259SPeter Brune #undef __FUNCT__
363c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown"
364c79ef259SPeter Brune /*@
365c79ef259SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
366c79ef259SPeter Brune    use on all levels.
367c79ef259SPeter Brune 
368c79ef259SPeter Brune    Logically Collective on PC
369c79ef259SPeter Brune 
370c79ef259SPeter Brune    Input Parameters:
371c79ef259SPeter Brune +  snes - the multigrid context
372c79ef259SPeter Brune -  n    - the number of smoothing steps
373c79ef259SPeter Brune 
374c79ef259SPeter Brune    Options Database Key:
375c79ef259SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
376c79ef259SPeter Brune 
377c79ef259SPeter Brune    Level: advanced
378c79ef259SPeter Brune 
379c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
380c79ef259SPeter Brune 
381c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
382c79ef259SPeter Brune @*/
383c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
384c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
385c79ef259SPeter Brune   PetscErrorCode ierr = 0;
386c79ef259SPeter Brune   PetscFunctionBegin;
387c79ef259SPeter Brune   fas->max_down_it = n;
388c79ef259SPeter Brune   if (fas->next) {
389c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
390c79ef259SPeter Brune   }
391c79ef259SPeter Brune   PetscFunctionReturn(0);
392c79ef259SPeter Brune }
393c79ef259SPeter Brune 
394c79ef259SPeter Brune #undef __FUNCT__
395421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
396c79ef259SPeter Brune /*@
397c79ef259SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
398c79ef259SPeter Brune    interpolation from l-1 to the lth level
399c79ef259SPeter Brune 
400c79ef259SPeter Brune    Logically Collective on PC and Mat
401c79ef259SPeter Brune 
402c79ef259SPeter Brune    Input Parameters:
403c79ef259SPeter Brune +  snes      - the multigrid context
404c79ef259SPeter Brune .  mat       - the interpolation operator
405c79ef259SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
406c79ef259SPeter Brune 
407c79ef259SPeter Brune    Level: advanced
408c79ef259SPeter Brune 
409c79ef259SPeter Brune    Notes:
410c79ef259SPeter Brune           Usually this is the same matrix used also to set the restriction
411c79ef259SPeter Brune     for the same level.
412c79ef259SPeter Brune 
413c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
414c79ef259SPeter Brune     out from the matrix size which one it is.
415c79ef259SPeter Brune 
416c79ef259SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
417c79ef259SPeter Brune 
418c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale()
419c79ef259SPeter Brune @*/
420421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
421421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
422d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
423421d9b32SPeter Brune 
424421d9b32SPeter Brune   PetscFunctionBegin;
425421d9b32SPeter Brune   if (level > top_level)
426421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
427421d9b32SPeter Brune   /* get to the correct level */
428d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
429421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
430421d9b32SPeter Brune   }
431421d9b32SPeter Brune   if (fas->level != level)
432421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
433421d9b32SPeter Brune   fas->interpolate = mat;
434421d9b32SPeter Brune   PetscFunctionReturn(0);
435421d9b32SPeter Brune }
436421d9b32SPeter Brune 
437421d9b32SPeter Brune #undef __FUNCT__
438421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
439c79ef259SPeter Brune /*@
440c79ef259SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
441c79ef259SPeter Brune    from level l to l-1.
442c79ef259SPeter Brune 
443c79ef259SPeter Brune    Logically Collective on SNES and Mat
444c79ef259SPeter Brune 
445c79ef259SPeter Brune    Input Parameters:
446c79ef259SPeter Brune +  snes  - the multigrid context
447c79ef259SPeter Brune .  mat   - the restriction matrix
448c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
449c79ef259SPeter Brune 
450c79ef259SPeter Brune    Level: advanced
451c79ef259SPeter Brune 
452c79ef259SPeter Brune    Notes:
453c79ef259SPeter Brune           Usually this is the same matrix used also to set the interpolation
454c79ef259SPeter Brune     for the same level.
455c79ef259SPeter Brune 
456c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
457c79ef259SPeter Brune     out from the matrix size which one it is.
458c79ef259SPeter Brune 
459c79ef259SPeter Brune          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
460c79ef259SPeter Brune     is used.
461c79ef259SPeter Brune 
462c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
463c79ef259SPeter Brune 
464c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
465c79ef259SPeter Brune @*/
466421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
467421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
468d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
469421d9b32SPeter Brune 
470421d9b32SPeter Brune   PetscFunctionBegin;
471421d9b32SPeter Brune   if (level > top_level)
472421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
473421d9b32SPeter Brune   /* get to the correct level */
474d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
475421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
476421d9b32SPeter Brune   }
477421d9b32SPeter Brune   if (fas->level != level)
478421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
479421d9b32SPeter Brune   fas->restrct = mat;
480421d9b32SPeter Brune   PetscFunctionReturn(0);
481421d9b32SPeter Brune }
482421d9b32SPeter Brune 
483efe1f98aSPeter Brune 
484421d9b32SPeter Brune #undef __FUNCT__
485421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
486c79ef259SPeter Brune /*@
487c79ef259SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
488c79ef259SPeter Brune    operator from level l to l-1.
489c79ef259SPeter Brune 
490c79ef259SPeter Brune    Logically Collective on SNES and Mat
491c79ef259SPeter Brune 
492c79ef259SPeter Brune    Input Parameters:
493c79ef259SPeter Brune +  snes   - the multigrid context
494c79ef259SPeter Brune .  rscale - the restriction scaling
495c79ef259SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
496c79ef259SPeter Brune 
497c79ef259SPeter Brune    Level: advanced
498c79ef259SPeter Brune 
499c79ef259SPeter Brune    Notes:
500c79ef259SPeter Brune          This is only used in the case that the injection is not set.
501c79ef259SPeter Brune 
502c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
503c79ef259SPeter Brune 
504c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
505c79ef259SPeter Brune @*/
506421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
507421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
508d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
509421d9b32SPeter Brune 
510421d9b32SPeter Brune   PetscFunctionBegin;
511421d9b32SPeter Brune   if (level > top_level)
512421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
513421d9b32SPeter Brune   /* get to the correct level */
514d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
515421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
516421d9b32SPeter Brune   }
517421d9b32SPeter Brune   if (fas->level != level)
518421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
519421d9b32SPeter Brune   fas->rscale = rscale;
520421d9b32SPeter Brune   PetscFunctionReturn(0);
521421d9b32SPeter Brune }
522421d9b32SPeter Brune 
523efe1f98aSPeter Brune 
524efe1f98aSPeter Brune #undef __FUNCT__
525efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection"
526c79ef259SPeter Brune /*@
527c79ef259SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
528c79ef259SPeter Brune    from level l to l-1.
529c79ef259SPeter Brune 
530c79ef259SPeter Brune    Logically Collective on SNES and Mat
531c79ef259SPeter Brune 
532c79ef259SPeter Brune    Input Parameters:
533c79ef259SPeter Brune +  snes  - the multigrid context
534c79ef259SPeter Brune .  mat   - the restriction matrix
535c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
536c79ef259SPeter Brune 
537c79ef259SPeter Brune    Level: advanced
538c79ef259SPeter Brune 
539c79ef259SPeter Brune    Notes:
540c79ef259SPeter Brune          If you do not set this, the restriction and rscale is used to
541c79ef259SPeter Brune    project the solution instead.
542c79ef259SPeter Brune 
543c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
544c79ef259SPeter Brune 
545c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
546c79ef259SPeter Brune @*/
547efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
548efe1f98aSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
549efe1f98aSPeter Brune   PetscInt top_level = fas->level,i;
550efe1f98aSPeter Brune 
551efe1f98aSPeter Brune   PetscFunctionBegin;
552efe1f98aSPeter Brune   if (level > top_level)
553efe1f98aSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
554efe1f98aSPeter Brune   /* get to the correct level */
555efe1f98aSPeter Brune   for (i = fas->level; i > level; i--) {
556efe1f98aSPeter Brune     fas = (SNES_FAS *)fas->next->data;
557efe1f98aSPeter Brune   }
558efe1f98aSPeter Brune   if (fas->level != level)
559efe1f98aSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
560efe1f98aSPeter Brune   fas->inject = mat;
561efe1f98aSPeter Brune   PetscFunctionReturn(0);
562efe1f98aSPeter Brune }
563efe1f98aSPeter Brune 
564421d9b32SPeter Brune #undef __FUNCT__
565421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
566421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
567421d9b32SPeter Brune {
56877df8cc4SPeter Brune   PetscErrorCode ierr = 0;
569421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
570421d9b32SPeter Brune 
571421d9b32SPeter Brune   PetscFunctionBegin;
572421d9b32SPeter Brune   /* destroy local data created in SNESSetup_FAS */
57351e86f29SPeter Brune #if 0
574293a7e31SPeter Brune   /* recurse -- reset should destroy the structures -- destroy should destroy the structures recursively */
57551e86f29SPeter Brune #endif
576ee1fd11aSPeter Brune   if (fas->next) ierr = SNESReset_FAS(fas->next);CHKERRQ(ierr);
57751e86f29SPeter Brune #if 0
57851e86f29SPeter Brune #endif
579421d9b32SPeter Brune   PetscFunctionReturn(0);
580421d9b32SPeter Brune }
581421d9b32SPeter Brune 
582421d9b32SPeter Brune #undef __FUNCT__
583421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
584421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
585421d9b32SPeter Brune {
586421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
587421d9b32SPeter Brune   PetscErrorCode ierr;
588421d9b32SPeter Brune 
589421d9b32SPeter Brune   PetscFunctionBegin;
590421d9b32SPeter Brune   /* recursively resets and then destroys */
591421d9b32SPeter Brune   ierr = SNESReset_FAS(snes);CHKERRQ(ierr);
59251e86f29SPeter Brune   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
59351e86f29SPeter Brune   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
594efe1f98aSPeter Brune   if (fas->inject) {
595efe1f98aSPeter Brune     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
596efe1f98aSPeter Brune   }
59751e86f29SPeter Brune   if (fas->interpolate == fas->restrct) {
59851e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
59951e86f29SPeter Brune     fas->restrct = PETSC_NULL;
60051e86f29SPeter Brune   } else {
60151e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
60251e86f29SPeter Brune     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
60351e86f29SPeter Brune   }
60451e86f29SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
60551e86f29SPeter Brune   if (snes->work)        ierr = VecDestroyVecs(snes->nwork,&snes->work);CHKERRQ(ierr);
606421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
607421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
608ee1fd11aSPeter Brune 
609421d9b32SPeter Brune   PetscFunctionReturn(0);
610421d9b32SPeter Brune }
611421d9b32SPeter Brune 
612421d9b32SPeter Brune #undef __FUNCT__
613421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
614421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
615421d9b32SPeter Brune {
6161a266240SBarry Smith   SNES_FAS       *fas = (SNES_FAS *) snes->data,*tmp;
617421d9b32SPeter Brune   PetscErrorCode ierr;
618efe1f98aSPeter Brune   VecScatter     injscatter;
619d1adcc6fSPeter Brune   PetscInt       dm_levels;
620d1adcc6fSPeter Brune 
621421d9b32SPeter Brune 
622421d9b32SPeter Brune   PetscFunctionBegin;
623eff52c0eSPeter Brune 
624cc05f883SPeter Brune   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
625d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
626d1adcc6fSPeter Brune     dm_levels++;
627cc05f883SPeter Brune     if (dm_levels > fas->levels) {
628d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
629cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
630d1adcc6fSPeter Brune     }
631d1adcc6fSPeter Brune   }
632d1adcc6fSPeter Brune 
633cc05f883SPeter Brune   if (!snes->work || snes->nwork != 2) {ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr);} /* work vectors used for intergrid transfers */
634cc05f883SPeter Brune 
635d1adcc6fSPeter Brune   /* should call the SNESSetFromOptions() only when appropriate */
6361a266240SBarry Smith   tmp = fas;
6371a266240SBarry Smith   while (tmp) {
638*6bed07a3SBarry Smith     if (tmp->upsmooth) {ierr = SNESSetFromOptions(tmp->upsmooth);CHKERRQ(ierr);}
639*6bed07a3SBarry Smith     if (tmp->downsmooth) {ierr = SNESSetFromOptions(tmp->downsmooth);CHKERRQ(ierr);}
6401a266240SBarry Smith     tmp = tmp->next ? (SNES_FAS*) tmp->next->data : 0;
6411a266240SBarry Smith   }
6421a266240SBarry Smith 
643421d9b32SPeter Brune   /* gets the solver ready for solution */
644646217ecSPeter Brune   if (fas->next) {
645646217ecSPeter Brune     if (snes->ops->computegs) {
646646217ecSPeter Brune       ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
647646217ecSPeter Brune     }
648646217ecSPeter Brune   }
649421d9b32SPeter Brune   if (snes->dm) {
650421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
651421d9b32SPeter Brune     if (fas->next) {
652421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
653ee78dd50SPeter Brune       if (!fas->next->dm) {
654ee78dd50SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
655ee78dd50SPeter Brune         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
656ee78dd50SPeter Brune       }
657421d9b32SPeter Brune       /* set the interpolation and restriction from the DM */
658ee78dd50SPeter Brune       if (!fas->interpolate) {
659421d9b32SPeter Brune         ierr = DMGetInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
660421d9b32SPeter Brune         fas->restrct = fas->interpolate;
661421d9b32SPeter Brune       }
662efe1f98aSPeter Brune       /* set the injection from the DM */
663efe1f98aSPeter Brune       if (!fas->inject) {
664efe1f98aSPeter Brune         ierr = DMGetInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
665efe1f98aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
666efe1f98aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
667efe1f98aSPeter Brune       }
668421d9b32SPeter Brune     }
669ee78dd50SPeter Brune     /* set the DMs of the pre and post-smoothers here */
670*6bed07a3SBarry Smith     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
671*6bed07a3SBarry Smith     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
672421d9b32SPeter Brune   }
673cc05f883SPeter Brune 
674fa9694d7SPeter Brune   if (fas->next) {
675fa9694d7SPeter Brune    /* gotta set up the solution vector for this to work */
676*6bed07a3SBarry Smith     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
677fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
678fa9694d7SPeter Brune   }
679fa9694d7SPeter Brune   /* got to set them all up at once */
680421d9b32SPeter Brune   PetscFunctionReturn(0);
681421d9b32SPeter Brune }
682421d9b32SPeter Brune 
683421d9b32SPeter Brune #undef __FUNCT__
684421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
685421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
686421d9b32SPeter Brune {
687ee78dd50SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
688ee78dd50SPeter Brune   PetscInt levels = 1;
689d1adcc6fSPeter Brune   PetscBool flg, smoothflg, smoothupflg, smoothdownflg, monflg;
690421d9b32SPeter Brune   PetscErrorCode ierr;
691ee78dd50SPeter Brune   const char * def_smooth = SNESNRICHARDSON;
692ee78dd50SPeter Brune   char pre_type[256];
693ee78dd50SPeter Brune   char post_type[256];
694ee78dd50SPeter Brune   char                    monfilename[PETSC_MAX_PATH_LEN];
695421d9b32SPeter Brune 
696421d9b32SPeter Brune   PetscFunctionBegin;
697c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
698ee78dd50SPeter Brune 
699ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
700ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
701ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
702c732cbdbSBarry Smith     if (!flg && snes->dm) {
703c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
704c732cbdbSBarry Smith       levels++;
705d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
706c732cbdbSBarry Smith     }
707ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
708ee78dd50SPeter Brune   }
709ee78dd50SPeter Brune 
710ee78dd50SPeter Brune   /* type of pre and/or post smoothers -- set both at once */
711ee78dd50SPeter Brune   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
712ee78dd50SPeter Brune   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
713d1adcc6fSPeter Brune   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
714d1adcc6fSPeter Brune   if (smoothflg) {
715ee78dd50SPeter Brune     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
716ee78dd50SPeter Brune   } else {
717d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr);
718d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr);
719ee78dd50SPeter Brune   }
720ee78dd50SPeter Brune 
721ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
722d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
723d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
724d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
725ee78dd50SPeter Brune 
726646217ecSPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
727ee78dd50SPeter Brune 
728eff52c0eSPeter Brune   ierr = PetscOptionsBool("-snes_fas_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr);
729eff52c0eSPeter Brune 
730ee78dd50SPeter Brune   /* other options for the coarsest level */
731ee78dd50SPeter Brune   if (fas->level == 0) {
732d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
733eff52c0eSPeter Brune     ierr = PetscOptionsBool("-snes_fas_coarse_ngs", "Use Nonlinear Gauss-Seidel smoother if provided", "SNESSetGS", fas->useGS, &fas->useGS, &flg);CHKERRQ(ierr);
734ee78dd50SPeter Brune   }
735ee78dd50SPeter Brune 
736421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
7378cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
738eff52c0eSPeter Brune 
739eff52c0eSPeter Brune   if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !fas->useGS)) {
7408cc86e31SPeter Brune     const char     *prefix;
7418cc86e31SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
7428cc86e31SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
7438cc86e31SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
7448cc86e31SPeter Brune     if (fas->level || (fas->levels == 1)) {
745eff52c0eSPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr);
7468cc86e31SPeter Brune     } else {
7478cc86e31SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
7488cc86e31SPeter Brune     }
7498cc86e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
7508cc86e31SPeter Brune     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
7518cc86e31SPeter Brune     if (snes->ops->computefunction) {
7528cc86e31SPeter Brune       ierr = SNESSetFunction(fas->downsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
7538cc86e31SPeter Brune     }
7548cc86e31SPeter Brune   }
7558cc86e31SPeter Brune 
756eff52c0eSPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !fas->useGS)) {
75767339d5cSBarry Smith     const char     *prefix;
75867339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
759ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
76067339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
761eff52c0eSPeter Brune     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr);
762293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
763ee78dd50SPeter Brune     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
76469bed9afSBarry Smith     if (snes->ops->computefunction) {
76569bed9afSBarry Smith       ierr = SNESSetFunction(fas->upsmooth,PETSC_NULL,snes->ops->computefunction,snes->funP);CHKERRQ(ierr);
76669bed9afSBarry Smith     }
767ee78dd50SPeter Brune   }
7681a266240SBarry Smith   if (fas->upsmooth) {
7691a266240SBarry Smith     ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
7701a266240SBarry Smith   }
7711a266240SBarry Smith 
7721a266240SBarry Smith   if (fas->downsmooth) {
7731a266240SBarry Smith     ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
7741a266240SBarry Smith   }
775ee78dd50SPeter Brune 
776ee78dd50SPeter Brune   if (monflg) {
777646217ecSPeter Brune     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
778794bee33SPeter Brune     /* set the monitors for the upsmoother and downsmoother */
779293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
780293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
781ee78dd50SPeter Brune   }
782ee78dd50SPeter Brune 
783ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
7841a266240SBarry Smith   if (fas->next) {ierr = SNESSetFromOptions_FAS(fas->next);CHKERRQ(ierr);}
785421d9b32SPeter Brune   PetscFunctionReturn(0);
786421d9b32SPeter Brune }
787421d9b32SPeter Brune 
788421d9b32SPeter Brune #undef __FUNCT__
789421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
790421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
791421d9b32SPeter Brune {
792421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
793421d9b32SPeter Brune   PetscBool      iascii;
794421d9b32SPeter Brune   PetscErrorCode ierr;
795421d9b32SPeter Brune   PetscInt levels = fas->levels;
796421d9b32SPeter Brune   PetscInt i;
797421d9b32SPeter Brune 
798421d9b32SPeter Brune   PetscFunctionBegin;
799421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
800421d9b32SPeter Brune   if (iascii) {
801421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
802421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
803421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
804421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
805ee78dd50SPeter Brune       if (fas->upsmooth) {
806ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
807421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
808ee78dd50SPeter Brune         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
809421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
810421d9b32SPeter Brune       } else {
811ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
812421d9b32SPeter Brune       }
813ee78dd50SPeter Brune       if (fas->downsmooth) {
814ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
815421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
816ee78dd50SPeter Brune         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
817421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
818421d9b32SPeter Brune       } else {
819ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
820421d9b32SPeter Brune       }
821eff52c0eSPeter Brune       if (fas->useGS) {
822eff52c0eSPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n",  fas->level);CHKERRQ(ierr);
823eff52c0eSPeter Brune       }
824421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
825421d9b32SPeter Brune     }
826421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
827421d9b32SPeter Brune   } else {
828421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
829421d9b32SPeter Brune   }
830421d9b32SPeter Brune   PetscFunctionReturn(0);
831421d9b32SPeter Brune }
832421d9b32SPeter Brune 
833421d9b32SPeter Brune /*
834fa9694d7SPeter Brune 
835fa9694d7SPeter Brune Defines the FAS cycle as:
836fa9694d7SPeter Brune 
837fa9694d7SPeter Brune fine problem: F(x) = 0
838fa9694d7SPeter Brune coarse problem: F^c(x) = b^c
839fa9694d7SPeter Brune 
840fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
841fa9694d7SPeter Brune 
842fa9694d7SPeter Brune correction:
843fa9694d7SPeter Brune 
844fa9694d7SPeter Brune x = x + I(x^c - Rx)
845fa9694d7SPeter Brune 
846421d9b32SPeter Brune  */
847fa9694d7SPeter Brune 
848421d9b32SPeter Brune #undef __FUNCT__
849421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private"
850646217ecSPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec X) {
851fa9694d7SPeter Brune 
852fa9694d7SPeter Brune   PetscErrorCode ierr;
853646217ecSPeter Brune   Vec X_c, Xo_c, F_c, B_c,F,B;
854fa9694d7SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
855646217ecSPeter Brune   PetscInt i, k;
8562d15c84fSPeter Brune   PetscReal fnorm;
857421d9b32SPeter Brune 
858421d9b32SPeter Brune   PetscFunctionBegin;
859f5a6d4f9SBarry Smith   F = snes->vec_func;
860646217ecSPeter Brune   B = snes->vec_rhs;
861fa9694d7SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
862d1adcc6fSPeter Brune   if (fas->downsmooth) {
863d1adcc6fSPeter Brune     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
864d1adcc6fSPeter Brune   } else if (snes->ops->computegs) {
865794bee33SPeter Brune     if (fas->monitor) {
866794bee33SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
867794bee33SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
868d1adcc6fSPeter Brune       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
869eff52c0eSPeter Brune       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
870d1adcc6fSPeter Brune       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
871794bee33SPeter Brune     }
872646217ecSPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
873646217ecSPeter Brune       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
874cc05f883SPeter Brune       if (fas->monitor) {
875794bee33SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
876794bee33SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
877d1adcc6fSPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
878eff52c0eSPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
879d1adcc6fSPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
880794bee33SPeter Brune       }
881646217ecSPeter Brune     }
882c90fad12SPeter Brune   } else if (snes->pc) {
883c90fad12SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
884fe6f9142SPeter Brune   }
885fa9694d7SPeter Brune   if (fas->next) {
886293a7e31SPeter Brune     for (i = 0; i < fas->n_cycles; i++) {
887ee78dd50SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
888794bee33SPeter Brune 
889c90fad12SPeter Brune       X_c  = fas->next->vec_sol;
890293a7e31SPeter Brune       Xo_c = fas->next->work[0];
891c90fad12SPeter Brune       F_c  = fas->next->vec_func;
892293a7e31SPeter Brune       B_c  = fas->next->work[1];
893efe1f98aSPeter Brune 
894efe1f98aSPeter Brune       /* inject the solution */
895efe1f98aSPeter Brune       if (fas->inject) {
896a3cb90a9SPeter Brune         ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
897efe1f98aSPeter Brune       } else {
898a3cb90a9SPeter Brune         ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
899a3cb90a9SPeter Brune         ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
900efe1f98aSPeter Brune       }
901293a7e31SPeter Brune       ierr = VecScale(F, -1.0);CHKERRQ(ierr);
902293a7e31SPeter Brune 
903293a7e31SPeter Brune       /* restrict the defect */
904293a7e31SPeter Brune       ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
905293a7e31SPeter Brune 
906c90fad12SPeter Brune       /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
907ee78dd50SPeter Brune       fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
908c90fad12SPeter Brune       ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
909293a7e31SPeter Brune       ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
910c90fad12SPeter Brune 
911ee78dd50SPeter Brune       /* set initial guess of the coarse problem to the projected fine solution */
912ee78dd50SPeter Brune       ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
913c90fad12SPeter Brune 
914c90fad12SPeter Brune       /* recurse to the next level */
915f5a6d4f9SBarry Smith       fas->next->vec_rhs = B_c;
916646217ecSPeter Brune       ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr);
917f5a6d4f9SBarry Smith       fas->next->vec_rhs = PETSC_NULL;
918ee78dd50SPeter Brune 
919fa9694d7SPeter Brune       /* correct as x <- x + I(x^c - Rx)*/
920fa9694d7SPeter Brune       ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
921ee78dd50SPeter Brune       ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr);
922293a7e31SPeter Brune     }
923fa9694d7SPeter Brune   }
924d1adcc6fSPeter Brune     /* up-smooth -- just update using the down-smoother */
925c90fad12SPeter Brune   if (fas->level != 0) {
926d1adcc6fSPeter Brune     if (fas->upsmooth) {
927d1adcc6fSPeter Brune       ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
928d1adcc6fSPeter Brune     } else if (snes->ops->computegs) {
929794bee33SPeter Brune       if (fas->monitor) {
930794bee33SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
931794bee33SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
932d1adcc6fSPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
933eff52c0eSPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
934d1adcc6fSPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
935794bee33SPeter Brune       }
936646217ecSPeter Brune       for (k = 0; k < fas->max_down_it; k++) {
937646217ecSPeter Brune         ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
938794bee33SPeter Brune         if (fas->monitor) {
939794bee33SPeter Brune           ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
940794bee33SPeter Brune           ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
941d1adcc6fSPeter Brune           ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
942eff52c0eSPeter Brune           ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
943d1adcc6fSPeter Brune           ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
944794bee33SPeter Brune         }
945646217ecSPeter Brune       }
946c90fad12SPeter Brune     } else if (snes->pc) {
947c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
948c90fad12SPeter Brune     }
949fe6f9142SPeter Brune   }
950fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
951fa9694d7SPeter Brune 
952fa9694d7SPeter Brune   PetscFunctionReturn(0);
953421d9b32SPeter Brune }
954421d9b32SPeter Brune 
955421d9b32SPeter Brune #undef __FUNCT__
956421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
957421d9b32SPeter Brune 
958421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
959421d9b32SPeter Brune {
960fa9694d7SPeter Brune   PetscErrorCode ierr;
961fe6f9142SPeter Brune   PetscInt i, maxits;
962f5a6d4f9SBarry Smith   Vec X, B,F;
963fe6f9142SPeter Brune   PetscReal fnorm;
964421d9b32SPeter Brune   PetscFunctionBegin;
965fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
966fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
967fa9694d7SPeter Brune   X = snes->vec_sol;
968fe6f9142SPeter Brune   B = snes->vec_rhs;
969f5a6d4f9SBarry Smith   F = snes->vec_func;
970293a7e31SPeter Brune 
971293a7e31SPeter Brune   /*norm setup */
972fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
973fe6f9142SPeter Brune   snes->iter = 0;
974fe6f9142SPeter Brune   snes->norm = 0.;
975fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
976fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
977fe6f9142SPeter Brune   if (snes->domainerror) {
978fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
979fe6f9142SPeter Brune     PetscFunctionReturn(0);
980fe6f9142SPeter Brune   }
981fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
982fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
983fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
984fe6f9142SPeter Brune   snes->norm = fnorm;
985fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
986fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
987fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
988fe6f9142SPeter Brune 
989fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
990fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
991fe6f9142SPeter Brune   /* test convergence */
992fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
993fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
994fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
995fe6f9142SPeter Brune     /* Call general purpose update function */
996646217ecSPeter Brune 
997fe6f9142SPeter Brune     if (snes->ops->update) {
998fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
999fe6f9142SPeter Brune     }
1000646217ecSPeter Brune     ierr = FASCycle_Private(snes, X);CHKERRQ(ierr);
1001c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1002c90fad12SPeter Brune     /* Monitor convergence */
1003c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1004c90fad12SPeter Brune     snes->iter = i+1;
1005c90fad12SPeter Brune     snes->norm = fnorm;
1006c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1007c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
1008c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1009c90fad12SPeter Brune     /* Test for convergence */
1010c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1011c90fad12SPeter Brune     if (snes->reason) break;
1012fe6f9142SPeter Brune   }
1013fe6f9142SPeter Brune   if (i == maxits) {
1014fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1015fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1016fe6f9142SPeter Brune   }
1017421d9b32SPeter Brune   PetscFunctionReturn(0);
1018421d9b32SPeter Brune }
1019