xref: /petsc/src/snes/impls/fas/fas.c (revision d28543b3ba94975fc020222886c1ad1afb4c04af)
1421d9b32SPeter Brune /* Defines the basic SNES object */
26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3421d9b32SPeter Brune 
4421d9b32SPeter Brune /*MC
5421d9b32SPeter Brune Full Approximation Scheme nonlinear multigrid solver.
6421d9b32SPeter Brune 
7421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
8421d9b32SPeter Brune 
9421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
10421d9b32SPeter Brune M*/
11421d9b32SPeter Brune 
12421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
13421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
14421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
15421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
16421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
17421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
18421d9b32SPeter Brune 
19421d9b32SPeter Brune EXTERN_C_BEGIN
20421d9b32SPeter Brune 
21421d9b32SPeter Brune #undef __FUNCT__
22421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
23421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes)
24421d9b32SPeter Brune {
25421d9b32SPeter Brune   SNES_FAS * fas;
26421d9b32SPeter Brune   PetscErrorCode ierr;
27421d9b32SPeter Brune 
28421d9b32SPeter Brune   PetscFunctionBegin;
29421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
30421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
31421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
32421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
33421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
34421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
35421d9b32SPeter Brune 
36ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
37ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
38ed020824SBarry Smith 
39421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
40421d9b32SPeter Brune   snes->data                  = (void*) fas;
41421d9b32SPeter Brune   fas->level                  = 0;
42293a7e31SPeter Brune   fas->levels                 = 1;
43ee78dd50SPeter Brune   fas->n_cycles               = 1;
44ee78dd50SPeter Brune   fas->max_up_it              = 1;
45ee78dd50SPeter Brune   fas->max_down_it            = 1;
46ee78dd50SPeter Brune   fas->upsmooth               = PETSC_NULL;
47ee78dd50SPeter Brune   fas->downsmooth             = PETSC_NULL;
48421d9b32SPeter Brune   fas->next                   = PETSC_NULL;
49421d9b32SPeter Brune   fas->interpolate            = PETSC_NULL;
50421d9b32SPeter Brune   fas->restrct                = PETSC_NULL;
51efe1f98aSPeter Brune   fas->inject                 = PETSC_NULL;
52cc05f883SPeter Brune   fas->monitor                = PETSC_NULL;
53cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
54421d9b32SPeter Brune   PetscFunctionReturn(0);
55421d9b32SPeter Brune }
56421d9b32SPeter Brune EXTERN_C_END
57421d9b32SPeter Brune 
58421d9b32SPeter Brune #undef __FUNCT__
59421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
60c79ef259SPeter Brune /*@
61722262beSPeter Brune    SNESFASGetLevels - Gets the number of levels in a FAS.
62c79ef259SPeter Brune 
63c79ef259SPeter Brune    Input Parameter:
64c79ef259SPeter Brune .  snes - the preconditioner context
65c79ef259SPeter Brune 
66c79ef259SPeter Brune    Output parameter:
67c79ef259SPeter Brune .  levels - the number of levels
68c79ef259SPeter Brune 
69c79ef259SPeter Brune    Level: advanced
70c79ef259SPeter Brune 
71c79ef259SPeter Brune .keywords: MG, get, levels, multigrid
72c79ef259SPeter Brune 
73c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels()
74c79ef259SPeter Brune @*/
75421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
76421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
77421d9b32SPeter Brune   PetscFunctionBegin;
78ee1fd11aSPeter Brune   *levels = fas->levels;
79421d9b32SPeter Brune   PetscFunctionReturn(0);
80421d9b32SPeter Brune }
81421d9b32SPeter Brune 
82421d9b32SPeter Brune #undef __FUNCT__
83646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles"
84c79ef259SPeter Brune /*@
85c79ef259SPeter Brune    SNESFASSetCycles - Sets the type cycles to use.  Use SNESFASSetCyclesOnLevel() for more
86c79ef259SPeter Brune    complicated cycling.
87c79ef259SPeter Brune 
88c79ef259SPeter Brune    Logically Collective on SNES
89c79ef259SPeter Brune 
90c79ef259SPeter Brune    Input Parameters:
91c79ef259SPeter Brune +  snes   - the multigrid context
92c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
93c79ef259SPeter Brune 
94c79ef259SPeter Brune    Options Database Key:
95c79ef259SPeter Brune $  -snes_fas_cycles 1 or 2
96c79ef259SPeter Brune 
97c79ef259SPeter Brune    Level: advanced
98c79ef259SPeter Brune 
99c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
100c79ef259SPeter Brune 
101c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel()
102c79ef259SPeter Brune @*/
103646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
104646217ecSPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
105646217ecSPeter Brune   PetscErrorCode ierr;
106646217ecSPeter Brune   PetscFunctionBegin;
107646217ecSPeter Brune   fas->n_cycles = cycles;
108646217ecSPeter Brune   if (fas->next) {
109646217ecSPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
110646217ecSPeter Brune   }
111646217ecSPeter Brune   PetscFunctionReturn(0);
112646217ecSPeter Brune }
113646217ecSPeter Brune 
114eff52c0eSPeter Brune #undef __FUNCT__
115c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel"
116c79ef259SPeter Brune /*@
117722262beSPeter Brune    SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level.
118c79ef259SPeter Brune 
119c79ef259SPeter Brune    Logically Collective on SNES
120c79ef259SPeter Brune 
121c79ef259SPeter Brune    Input Parameters:
122c79ef259SPeter Brune +  snes   - the multigrid context
123c79ef259SPeter Brune .  level  - the level to set the number of cycles on
124c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
125c79ef259SPeter Brune 
126c79ef259SPeter Brune    Level: advanced
127c79ef259SPeter Brune 
128c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
129c79ef259SPeter Brune 
130c79ef259SPeter Brune .seealso: SNESFASSetCycles()
131c79ef259SPeter Brune @*/
132c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) {
133c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
134c79ef259SPeter Brune   PetscInt top_level = fas->level,i;
135c79ef259SPeter Brune 
136c79ef259SPeter Brune   PetscFunctionBegin;
137c79ef259SPeter Brune   if (level > top_level)
138c79ef259SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
139c79ef259SPeter Brune   /* get to the correct level */
140c79ef259SPeter Brune   for (i = fas->level; i > level; i--) {
141c79ef259SPeter Brune     fas = (SNES_FAS *)fas->next->data;
142c79ef259SPeter Brune   }
143c79ef259SPeter Brune   if (fas->level != level)
144c79ef259SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
145c79ef259SPeter Brune   fas->n_cycles = cycles;
146c79ef259SPeter Brune   PetscFunctionReturn(0);
147c79ef259SPeter Brune }
148c79ef259SPeter Brune 
149c79ef259SPeter Brune #undef __FUNCT__
150eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS"
151aeed3662SMatthew G Knepley /*@C
152c79ef259SPeter Brune    SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used.
153c79ef259SPeter Brune    Use SNESFASSetGSOnLevel() for more complicated staging of smoothers
154c79ef259SPeter Brune    and nonlinear preconditioners.
155c79ef259SPeter Brune 
156c79ef259SPeter Brune    Logically Collective on SNES
157c79ef259SPeter Brune 
158c79ef259SPeter Brune    Input Parameters:
159c79ef259SPeter Brune +  snes    - the multigrid context
160c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
161c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
162c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
163c79ef259SPeter Brune 
164c79ef259SPeter Brune    Level: advanced
165c79ef259SPeter Brune 
166c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid
167c79ef259SPeter Brune 
168c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel()
169c79ef259SPeter Brune @*/
170eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
171eff52c0eSPeter Brune   PetscErrorCode ierr = 0;
172eff52c0eSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
173eff52c0eSPeter Brune   PetscFunctionBegin;
174eff52c0eSPeter Brune 
175eff52c0eSPeter Brune   /* use or don't use it according to user wishes*/
176d28543b3SPeter Brune   snes->usegs = use_gs;
177eff52c0eSPeter Brune   if (gsfunc) {
178eff52c0eSPeter Brune     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
179eff52c0eSPeter Brune     /* push the provided GS up the tree */
180eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
181eff52c0eSPeter Brune   } else if (snes->ops->computegs) {
182eff52c0eSPeter Brune     /* assume that the user has set the GS solver at this level */
183eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
184eff52c0eSPeter Brune   } else if (use_gs) {
185eff52c0eSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level);
186eff52c0eSPeter Brune   }
187eff52c0eSPeter Brune   PetscFunctionReturn(0);
188eff52c0eSPeter Brune }
189eff52c0eSPeter Brune 
190eff52c0eSPeter Brune #undef __FUNCT__
191eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel"
192aeed3662SMatthew G Knepley /*@C
193c79ef259SPeter Brune    SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level.
194c79ef259SPeter Brune 
195c79ef259SPeter Brune    Logically Collective on SNES
196c79ef259SPeter Brune 
197c79ef259SPeter Brune    Input Parameters:
198c79ef259SPeter Brune +  snes    - the multigrid context
199c79ef259SPeter Brune .  level   - the level to set the nonlinear smoother on
200c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
201c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
202c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
203c79ef259SPeter Brune 
204c79ef259SPeter Brune    Level: advanced
205c79ef259SPeter Brune 
206c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
207c79ef259SPeter Brune 
208c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS()
209c79ef259SPeter Brune @*/
210eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
211eff52c0eSPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
212eff52c0eSPeter Brune   PetscErrorCode ierr;
213eff52c0eSPeter Brune   PetscInt       top_level = fas->level,i;
214eff52c0eSPeter Brune   SNES           cur_snes = snes;
215eff52c0eSPeter Brune   PetscFunctionBegin;
216eff52c0eSPeter Brune   if (level > top_level)
217eff52c0eSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
218eff52c0eSPeter Brune   /* get to the correct level */
219eff52c0eSPeter Brune   for (i = fas->level; i > level; i--) {
220eff52c0eSPeter Brune     fas = (SNES_FAS *)fas->next->data;
221eff52c0eSPeter Brune     cur_snes = fas->next;
222eff52c0eSPeter Brune   }
223eff52c0eSPeter Brune   if (fas->level != level)
224eff52c0eSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
225d28543b3SPeter Brune   snes->usegs = use_gs;
226eff52c0eSPeter Brune   if (gsfunc) {
227eff52c0eSPeter Brune     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
228eff52c0eSPeter Brune   }
229eff52c0eSPeter Brune   PetscFunctionReturn(0);
230eff52c0eSPeter Brune }
231eff52c0eSPeter Brune 
232646217ecSPeter Brune #undef __FUNCT__
233421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
234c79ef259SPeter Brune /*@
235c79ef259SPeter Brune    SNESFASGetSNES - Gets the SNES corresponding to a particular
236c79ef259SPeter Brune    level of the FAS hierarchy.
237c79ef259SPeter Brune 
238c79ef259SPeter Brune    Input Parameters:
239c79ef259SPeter Brune +  snes    - the multigrid context
240c79ef259SPeter Brune    level   - the level to get
241c79ef259SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
242c79ef259SPeter Brune 
243c79ef259SPeter Brune    Level: advanced
244c79ef259SPeter Brune 
245c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
246c79ef259SPeter Brune 
247c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
248c79ef259SPeter Brune @*/
249421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
250421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
251421d9b32SPeter Brune   PetscInt levels = fas->level;
252421d9b32SPeter Brune   PetscInt i;
253421d9b32SPeter Brune   PetscFunctionBegin;
254421d9b32SPeter Brune   *lsnes = snes;
255421d9b32SPeter Brune   if (fas->level < level) {
256421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
257421d9b32SPeter Brune   }
258421d9b32SPeter Brune   if (level > levels - 1) {
259421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
260421d9b32SPeter Brune   }
261421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
262421d9b32SPeter Brune     *lsnes = fas->next;
263421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
264421d9b32SPeter Brune   }
265421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
266421d9b32SPeter Brune   PetscFunctionReturn(0);
267421d9b32SPeter Brune }
268421d9b32SPeter Brune 
269421d9b32SPeter Brune #undef __FUNCT__
270421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
271aeed3662SMatthew G Knepley /*@C
272c79ef259SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
273c79ef259SPeter Brune    Must be called before any other FAS routine.
274c79ef259SPeter Brune 
275c79ef259SPeter Brune    Input Parameters:
276c79ef259SPeter Brune +  snes   - the snes context
277c79ef259SPeter Brune .  levels - the number of levels
278c79ef259SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
279c79ef259SPeter Brune             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
280c79ef259SPeter Brune             Fortran.
281c79ef259SPeter Brune 
282c79ef259SPeter Brune    Level: intermediate
283c79ef259SPeter Brune 
284c79ef259SPeter Brune    Notes:
285c79ef259SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
286c79ef259SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
287c79ef259SPeter Brune 
288c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid
289c79ef259SPeter Brune 
290c79ef259SPeter Brune .seealso: SNESFASGetLevels()
291c79ef259SPeter Brune @*/
292421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
293421d9b32SPeter Brune   PetscErrorCode ierr;
294421d9b32SPeter Brune   PetscInt i;
295421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
296421d9b32SPeter Brune   MPI_Comm comm;
297421d9b32SPeter Brune   PetscFunctionBegin;
298ee1fd11aSPeter Brune    comm = ((PetscObject)snes)->comm;
299ee1fd11aSPeter Brune   if (levels == fas->levels) {
300ee1fd11aSPeter Brune     if (!comms)
301ee1fd11aSPeter Brune       PetscFunctionReturn(0);
302ee1fd11aSPeter Brune   }
303421d9b32SPeter Brune   /* user has changed the number of levels; reset */
304421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
305421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
306421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
307ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
308421d9b32SPeter Brune   /* setup the finest level */
309421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
310421d9b32SPeter Brune     if (comms) comm = comms[i];
311421d9b32SPeter Brune     fas->level = i;
312421d9b32SPeter Brune     fas->levels = levels;
313ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
314421d9b32SPeter Brune     if (i > 0) {
315421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
3164a6b3fb8SBarry Smith       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
317421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
318794bee33SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
319421d9b32SPeter Brune       fas = (SNES_FAS *)fas->next->data;
320421d9b32SPeter Brune     }
321421d9b32SPeter Brune   }
322421d9b32SPeter Brune   PetscFunctionReturn(0);
323421d9b32SPeter Brune }
324421d9b32SPeter Brune 
325421d9b32SPeter Brune #undef __FUNCT__
326c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp"
327c79ef259SPeter Brune /*@
328c79ef259SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
329c79ef259SPeter Brune    use on all levels.
330c79ef259SPeter Brune 
331c79ef259SPeter Brune    Logically Collective on PC
332c79ef259SPeter Brune 
333c79ef259SPeter Brune    Input Parameters:
334c79ef259SPeter Brune +  snes - the multigrid context
335c79ef259SPeter Brune -  n    - the number of smoothing steps
336c79ef259SPeter Brune 
337c79ef259SPeter Brune    Options Database Key:
338d28543b3SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
339c79ef259SPeter Brune 
340c79ef259SPeter Brune    Level: advanced
341c79ef259SPeter Brune 
342c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
343c79ef259SPeter Brune 
344c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
345c79ef259SPeter Brune @*/
346c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
347c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
348c79ef259SPeter Brune   PetscErrorCode ierr = 0;
349c79ef259SPeter Brune   PetscFunctionBegin;
350c79ef259SPeter Brune   fas->max_up_it = n;
351c79ef259SPeter Brune   if (fas->next) {
352c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
353c79ef259SPeter Brune   }
354c79ef259SPeter Brune   PetscFunctionReturn(0);
355c79ef259SPeter Brune }
356c79ef259SPeter Brune 
357c79ef259SPeter Brune #undef __FUNCT__
358c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown"
359c79ef259SPeter Brune /*@
360c79ef259SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
361c79ef259SPeter Brune    use on all levels.
362c79ef259SPeter Brune 
363c79ef259SPeter Brune    Logically Collective on PC
364c79ef259SPeter Brune 
365c79ef259SPeter Brune    Input Parameters:
366c79ef259SPeter Brune +  snes - the multigrid context
367c79ef259SPeter Brune -  n    - the number of smoothing steps
368c79ef259SPeter Brune 
369c79ef259SPeter Brune    Options Database Key:
370d28543b3SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
371c79ef259SPeter Brune 
372c79ef259SPeter Brune    Level: advanced
373c79ef259SPeter Brune 
374c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
375c79ef259SPeter Brune 
376c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
377c79ef259SPeter Brune @*/
378c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
379c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
380c79ef259SPeter Brune   PetscErrorCode ierr = 0;
381c79ef259SPeter Brune   PetscFunctionBegin;
382c79ef259SPeter Brune   fas->max_down_it = n;
383c79ef259SPeter Brune   if (fas->next) {
384c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
385c79ef259SPeter Brune   }
386c79ef259SPeter Brune   PetscFunctionReturn(0);
387c79ef259SPeter Brune }
388c79ef259SPeter Brune 
389c79ef259SPeter Brune #undef __FUNCT__
390421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
391c79ef259SPeter Brune /*@
392c79ef259SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
393c79ef259SPeter Brune    interpolation from l-1 to the lth level
394c79ef259SPeter Brune 
395c79ef259SPeter Brune    Input Parameters:
396c79ef259SPeter Brune +  snes      - the multigrid context
397c79ef259SPeter Brune .  mat       - the interpolation operator
398c79ef259SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
399c79ef259SPeter Brune 
400c79ef259SPeter Brune    Level: advanced
401c79ef259SPeter Brune 
402c79ef259SPeter Brune    Notes:
403c79ef259SPeter Brune           Usually this is the same matrix used also to set the restriction
404c79ef259SPeter Brune     for the same level.
405c79ef259SPeter Brune 
406c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
407c79ef259SPeter Brune     out from the matrix size which one it is.
408c79ef259SPeter Brune 
409c79ef259SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
410c79ef259SPeter Brune 
411c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale()
412c79ef259SPeter Brune @*/
413421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
414421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
415d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
416421d9b32SPeter Brune 
417421d9b32SPeter Brune   PetscFunctionBegin;
418421d9b32SPeter Brune   if (level > top_level)
419421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
420421d9b32SPeter Brune   /* get to the correct level */
421d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
422421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
423421d9b32SPeter Brune   }
424421d9b32SPeter Brune   if (fas->level != level)
425421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
426421d9b32SPeter Brune   fas->interpolate = mat;
427421d9b32SPeter Brune   PetscFunctionReturn(0);
428421d9b32SPeter Brune }
429421d9b32SPeter Brune 
430421d9b32SPeter Brune #undef __FUNCT__
431421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
432c79ef259SPeter Brune /*@
433c79ef259SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
434c79ef259SPeter Brune    from level l to l-1.
435c79ef259SPeter Brune 
436c79ef259SPeter Brune    Input Parameters:
437c79ef259SPeter Brune +  snes  - the multigrid context
438c79ef259SPeter Brune .  mat   - the restriction matrix
439c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
440c79ef259SPeter Brune 
441c79ef259SPeter Brune    Level: advanced
442c79ef259SPeter Brune 
443c79ef259SPeter Brune    Notes:
444c79ef259SPeter Brune           Usually this is the same matrix used also to set the interpolation
445c79ef259SPeter Brune     for the same level.
446c79ef259SPeter Brune 
447c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
448c79ef259SPeter Brune     out from the matrix size which one it is.
449c79ef259SPeter Brune 
450c79ef259SPeter Brune          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
451c79ef259SPeter Brune     is used.
452c79ef259SPeter Brune 
453c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
454c79ef259SPeter Brune 
455c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
456c79ef259SPeter Brune @*/
457421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
458421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
459d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
460421d9b32SPeter Brune 
461421d9b32SPeter Brune   PetscFunctionBegin;
462421d9b32SPeter Brune   if (level > top_level)
463421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
464421d9b32SPeter Brune   /* get to the correct level */
465d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
466421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
467421d9b32SPeter Brune   }
468421d9b32SPeter Brune   if (fas->level != level)
469421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
470421d9b32SPeter Brune   fas->restrct = mat;
471421d9b32SPeter Brune   PetscFunctionReturn(0);
472421d9b32SPeter Brune }
473421d9b32SPeter Brune 
474efe1f98aSPeter Brune 
475421d9b32SPeter Brune #undef __FUNCT__
476421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
477c79ef259SPeter Brune /*@
478c79ef259SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
479c79ef259SPeter Brune    operator from level l to l-1.
480c79ef259SPeter Brune 
481c79ef259SPeter Brune    Input Parameters:
482c79ef259SPeter Brune +  snes   - the multigrid context
483c79ef259SPeter Brune .  rscale - the restriction scaling
484c79ef259SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
485c79ef259SPeter Brune 
486c79ef259SPeter Brune    Level: advanced
487c79ef259SPeter Brune 
488c79ef259SPeter Brune    Notes:
489c79ef259SPeter Brune          This is only used in the case that the injection is not set.
490c79ef259SPeter Brune 
491c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
492c79ef259SPeter Brune 
493c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
494c79ef259SPeter Brune @*/
495421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
496421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
497d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
498421d9b32SPeter Brune 
499421d9b32SPeter Brune   PetscFunctionBegin;
500421d9b32SPeter Brune   if (level > top_level)
501421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
502421d9b32SPeter Brune   /* get to the correct level */
503d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
504421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
505421d9b32SPeter Brune   }
506421d9b32SPeter Brune   if (fas->level != level)
507421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
508421d9b32SPeter Brune   fas->rscale = rscale;
509421d9b32SPeter Brune   PetscFunctionReturn(0);
510421d9b32SPeter Brune }
511421d9b32SPeter Brune 
512efe1f98aSPeter Brune 
513efe1f98aSPeter Brune #undef __FUNCT__
514efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection"
515c79ef259SPeter Brune /*@
516c79ef259SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
517c79ef259SPeter Brune    from level l to l-1.
518c79ef259SPeter Brune 
519c79ef259SPeter Brune    Input Parameters:
520c79ef259SPeter Brune +  snes  - the multigrid context
521c79ef259SPeter Brune .  mat   - the restriction matrix
522c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
523c79ef259SPeter Brune 
524c79ef259SPeter Brune    Level: advanced
525c79ef259SPeter Brune 
526c79ef259SPeter Brune    Notes:
527c79ef259SPeter Brune          If you do not set this, the restriction and rscale is used to
528c79ef259SPeter Brune    project the solution instead.
529c79ef259SPeter Brune 
530c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
531c79ef259SPeter Brune 
532c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
533c79ef259SPeter Brune @*/
534efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
535efe1f98aSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
536efe1f98aSPeter Brune   PetscInt top_level = fas->level,i;
537efe1f98aSPeter Brune 
538efe1f98aSPeter Brune   PetscFunctionBegin;
539efe1f98aSPeter Brune   if (level > top_level)
540efe1f98aSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
541efe1f98aSPeter Brune   /* get to the correct level */
542efe1f98aSPeter Brune   for (i = fas->level; i > level; i--) {
543efe1f98aSPeter Brune     fas = (SNES_FAS *)fas->next->data;
544efe1f98aSPeter Brune   }
545efe1f98aSPeter Brune   if (fas->level != level)
546efe1f98aSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
547efe1f98aSPeter Brune   fas->inject = mat;
548efe1f98aSPeter Brune   PetscFunctionReturn(0);
549efe1f98aSPeter Brune }
550efe1f98aSPeter Brune 
551421d9b32SPeter Brune #undef __FUNCT__
552421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
553421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
554421d9b32SPeter Brune {
55577df8cc4SPeter Brune   PetscErrorCode ierr = 0;
556421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
557421d9b32SPeter Brune 
558421d9b32SPeter Brune   PetscFunctionBegin;
559742fe5e2SPeter Brune   if (fas->upsmooth)   ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
560742fe5e2SPeter Brune   if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
561742fe5e2SPeter Brune   if (fas->next)       ierr = SNESReset(fas->next);CHKERRQ(ierr);
562421d9b32SPeter Brune   PetscFunctionReturn(0);
563421d9b32SPeter Brune }
564421d9b32SPeter Brune 
565421d9b32SPeter Brune #undef __FUNCT__
566421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
567421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
568421d9b32SPeter Brune {
569421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
570742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
571421d9b32SPeter Brune 
572421d9b32SPeter Brune   PetscFunctionBegin;
573421d9b32SPeter Brune   /* recursively resets and then destroys */
57451e86f29SPeter Brune   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
57551e86f29SPeter Brune   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
576efe1f98aSPeter Brune   if (fas->inject) {
577efe1f98aSPeter Brune     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
578efe1f98aSPeter Brune   }
57951e86f29SPeter Brune   if (fas->interpolate == fas->restrct) {
58051e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
58151e86f29SPeter Brune     fas->restrct = PETSC_NULL;
58251e86f29SPeter Brune   } else {
58351e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
58451e86f29SPeter Brune     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
58551e86f29SPeter Brune   }
58651e86f29SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
587421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
588421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
589ee1fd11aSPeter Brune 
590421d9b32SPeter Brune   PetscFunctionReturn(0);
591421d9b32SPeter Brune }
592421d9b32SPeter Brune 
593421d9b32SPeter Brune #undef __FUNCT__
594421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
595421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
596421d9b32SPeter Brune {
59748bfdf8aSPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
598421d9b32SPeter Brune   PetscErrorCode ierr;
599efe1f98aSPeter Brune   VecScatter     injscatter;
600d1adcc6fSPeter Brune   PetscInt       dm_levels;
601d1adcc6fSPeter Brune 
602421d9b32SPeter Brune 
603421d9b32SPeter Brune   PetscFunctionBegin;
604eff52c0eSPeter Brune 
605cc05f883SPeter Brune   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
606d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
607d1adcc6fSPeter Brune     dm_levels++;
608cc05f883SPeter Brune     if (dm_levels > fas->levels) {
609d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
610cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
611d1adcc6fSPeter Brune     }
612d1adcc6fSPeter Brune   }
613d1adcc6fSPeter Brune 
6142f7ea302SPeter Brune   ierr = SNESDefaultGetWork(snes, 1);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
615cc05f883SPeter Brune 
61648bfdf8aSPeter Brune   /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */
61748bfdf8aSPeter Brune   if (fas->upsmooth) {
61848bfdf8aSPeter Brune     ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);
61948bfdf8aSPeter Brune     if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) {
62048bfdf8aSPeter Brune       ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
62148bfdf8aSPeter Brune     }
62248bfdf8aSPeter Brune     if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) {
62348bfdf8aSPeter Brune       ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
62448bfdf8aSPeter Brune     }
62548bfdf8aSPeter Brune    if (snes->ops->computegs && !fas->upsmooth->ops->computegs) {
62648bfdf8aSPeter Brune       ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
62748bfdf8aSPeter Brune     }
62848bfdf8aSPeter Brune   }
62948bfdf8aSPeter Brune   if (fas->downsmooth) {
63048bfdf8aSPeter Brune     ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);
63148bfdf8aSPeter Brune     if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) {
63248bfdf8aSPeter Brune       ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
63348bfdf8aSPeter Brune     }
63448bfdf8aSPeter Brune     if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) {
63548bfdf8aSPeter Brune       ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
63648bfdf8aSPeter Brune     }
63748bfdf8aSPeter Brune    if (snes->ops->computegs && !fas->downsmooth->ops->computegs) {
63848bfdf8aSPeter Brune       ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
63948bfdf8aSPeter Brune     }
6401a266240SBarry Smith   }
64148bfdf8aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
642646217ecSPeter Brune   if (fas->next) {
64348bfdf8aSPeter Brune     if (snes->ops->computefunction && !fas->next->ops->computefunction) {
64448bfdf8aSPeter Brune       ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
64548bfdf8aSPeter Brune     }
64648bfdf8aSPeter Brune     if (snes->ops->computejacobian && !fas->next->ops->computejacobian) {
64748bfdf8aSPeter Brune       ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
64848bfdf8aSPeter Brune     }
64948bfdf8aSPeter Brune    if (snes->ops->computegs && !fas->next->ops->computegs) {
650646217ecSPeter Brune       ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
651646217ecSPeter Brune     }
652646217ecSPeter Brune   }
653421d9b32SPeter Brune   if (snes->dm) {
654421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
655421d9b32SPeter Brune     if (fas->next) {
656421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
657ee78dd50SPeter Brune       if (!fas->next->dm) {
658ee78dd50SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
659ee78dd50SPeter Brune         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
660ee78dd50SPeter Brune       }
661421d9b32SPeter Brune       /* set the interpolation and restriction from the DM */
662ee78dd50SPeter Brune       if (!fas->interpolate) {
663e727c939SJed Brown         ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
664421d9b32SPeter Brune         fas->restrct = fas->interpolate;
665421d9b32SPeter Brune       }
666efe1f98aSPeter Brune       /* set the injection from the DM */
667efe1f98aSPeter Brune       if (!fas->inject) {
668e727c939SJed Brown         ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
669efe1f98aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
670efe1f98aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
671efe1f98aSPeter Brune       }
672421d9b32SPeter Brune     }
673ee78dd50SPeter Brune     /* set the DMs of the pre and post-smoothers here */
674*6bed07a3SBarry Smith     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
675*6bed07a3SBarry Smith     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
676421d9b32SPeter Brune   }
677cc05f883SPeter Brune 
678fa9694d7SPeter Brune   if (fas->next) {
679fa9694d7SPeter Brune    /* gotta set up the solution vector for this to work */
680*6bed07a3SBarry Smith     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
681742fe5e2SPeter Brune     if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);}
682fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
683fa9694d7SPeter Brune   }
684fa9694d7SPeter Brune   /* got to set them all up at once */
685421d9b32SPeter Brune   PetscFunctionReturn(0);
686421d9b32SPeter Brune }
687421d9b32SPeter Brune 
688421d9b32SPeter Brune #undef __FUNCT__
689421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
690421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
691421d9b32SPeter Brune {
692ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
693ee78dd50SPeter Brune   PetscInt       levels = 1;
694d1adcc6fSPeter Brune   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, monflg;
695421d9b32SPeter Brune   PetscErrorCode ierr;
696ee78dd50SPeter Brune   const char     *def_smooth = SNESNRICHARDSON;
697ee78dd50SPeter Brune   char           pre_type[256];
698ee78dd50SPeter Brune   char           post_type[256];
699ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
700421d9b32SPeter Brune 
701421d9b32SPeter Brune   PetscFunctionBegin;
702c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
703ee78dd50SPeter Brune 
704ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
705ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
706ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
707c732cbdbSBarry Smith     if (!flg && snes->dm) {
708c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
709c732cbdbSBarry Smith       levels++;
710d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
711c732cbdbSBarry Smith     }
712ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
713ee78dd50SPeter Brune   }
714ee78dd50SPeter Brune 
715ee78dd50SPeter Brune   /* type of pre and/or post smoothers -- set both at once */
716ee78dd50SPeter Brune   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
717ee78dd50SPeter Brune   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
718d1adcc6fSPeter Brune   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
719d1adcc6fSPeter Brune   if (smoothflg) {
720ee78dd50SPeter Brune     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
721ee78dd50SPeter Brune   } else {
722d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr);
723d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr);
724ee78dd50SPeter Brune   }
725ee78dd50SPeter Brune 
726ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
727d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
728d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
729d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
730ee78dd50SPeter Brune 
731646217ecSPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
732ee78dd50SPeter Brune 
733ee78dd50SPeter Brune   /* other options for the coarsest level */
734ee78dd50SPeter Brune   if (fas->level == 0) {
735d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
736ee78dd50SPeter Brune   }
737ee78dd50SPeter Brune 
738421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
7398cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
740eff52c0eSPeter Brune 
741d28543b3SPeter Brune   if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->usegs)) {
7428cc86e31SPeter Brune     const char     *prefix;
7438cc86e31SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
7448cc86e31SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
7458cc86e31SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
7468cc86e31SPeter Brune     if (fas->level || (fas->levels == 1)) {
747eff52c0eSPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr);
7488cc86e31SPeter Brune     } else {
7498cc86e31SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
7508cc86e31SPeter Brune     }
7518cc86e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
7528cc86e31SPeter Brune     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
7538cc86e31SPeter Brune   }
7548cc86e31SPeter Brune 
755d28543b3SPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->usegs)) {
75667339d5cSBarry Smith     const char     *prefix;
75767339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
758ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
75967339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
760eff52c0eSPeter Brune     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr);
761293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
762ee78dd50SPeter Brune     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
763ee78dd50SPeter Brune   }
7641a266240SBarry Smith   if (fas->upsmooth) {
7651a266240SBarry Smith     ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
7661a266240SBarry Smith   }
7671a266240SBarry Smith 
7681a266240SBarry Smith   if (fas->downsmooth) {
7691a266240SBarry Smith     ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
7701a266240SBarry Smith   }
771ee78dd50SPeter Brune 
772742fe5e2SPeter Brune   if (fas->level != fas->levels - 1) {
773742fe5e2SPeter Brune     ierr = SNESSetTolerances(snes, 0.0, 0.0, 0.0, fas->n_cycles, 1000);CHKERRQ(ierr);
774742fe5e2SPeter Brune   }
775742fe5e2SPeter 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 */
7792f7ea302SPeter Brune     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
780742fe5e2SPeter Brune     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
781293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
782293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
783d28543b3SPeter Brune   } else {
784d28543b3SPeter Brune     /* unset the monitors on the coarse levels */
785d28543b3SPeter Brune     if (fas->level != fas->levels - 1) {
786d28543b3SPeter Brune       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
787d28543b3SPeter Brune     }
788ee78dd50SPeter Brune   }
789ee78dd50SPeter Brune 
790ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
791d28543b3SPeter Brune   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
792421d9b32SPeter Brune   PetscFunctionReturn(0);
793421d9b32SPeter Brune }
794421d9b32SPeter Brune 
795421d9b32SPeter Brune #undef __FUNCT__
796421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
797421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
798421d9b32SPeter Brune {
799421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
800421d9b32SPeter Brune   PetscBool      iascii;
801421d9b32SPeter Brune   PetscErrorCode ierr;
802421d9b32SPeter Brune   PetscInt levels = fas->levels;
803421d9b32SPeter Brune   PetscInt i;
804421d9b32SPeter Brune 
805421d9b32SPeter Brune   PetscFunctionBegin;
806421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
807421d9b32SPeter Brune   if (iascii) {
808421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
809421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
810421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
811421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
812ee78dd50SPeter Brune       if (fas->upsmooth) {
813ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
814421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
815ee78dd50SPeter Brune         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
816421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
817421d9b32SPeter Brune       } else {
818ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
819421d9b32SPeter Brune       }
820ee78dd50SPeter Brune       if (fas->downsmooth) {
821ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
822421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
823ee78dd50SPeter Brune         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
824421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
825421d9b32SPeter Brune       } else {
826ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
827421d9b32SPeter Brune       }
828d28543b3SPeter Brune       if (snes->usegs) {
829eff52c0eSPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n",  fas->level);CHKERRQ(ierr);
830eff52c0eSPeter Brune       }
831421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
832421d9b32SPeter Brune     }
833421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
834421d9b32SPeter Brune   } else {
835421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
836421d9b32SPeter Brune   }
837421d9b32SPeter Brune   PetscFunctionReturn(0);
838421d9b32SPeter Brune }
839421d9b32SPeter Brune 
840421d9b32SPeter Brune /*
841fa9694d7SPeter Brune 
842fa9694d7SPeter Brune Defines the FAS cycle as:
843fa9694d7SPeter Brune 
844fa9694d7SPeter Brune fine problem: F(x) = 0
845fa9694d7SPeter Brune coarse problem: F^c(x) = b^c
846fa9694d7SPeter Brune 
847fa9694d7SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
848fa9694d7SPeter Brune 
849fa9694d7SPeter Brune correction:
850fa9694d7SPeter Brune 
851fa9694d7SPeter Brune x = x + I(x^c - Rx)
852fa9694d7SPeter Brune 
853421d9b32SPeter Brune  */
854fa9694d7SPeter Brune 
855421d9b32SPeter Brune #undef __FUNCT__
856421d9b32SPeter Brune #define __FUNCT__ "FASCycle_Private"
857646217ecSPeter Brune PetscErrorCode FASCycle_Private(SNES snes, Vec X) {
858fa9694d7SPeter Brune 
859fa9694d7SPeter Brune   PetscErrorCode      ierr;
860646217ecSPeter Brune   Vec                 X_c, Xo_c, F_c, B_c,F,B;
861fa9694d7SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
862742fe5e2SPeter Brune   PetscInt            k;
8632d15c84fSPeter Brune   PetscReal           fnorm;
864742fe5e2SPeter Brune   SNESConvergedReason reason;
865421d9b32SPeter Brune 
866421d9b32SPeter Brune   PetscFunctionBegin;
867f5a6d4f9SBarry Smith   F = snes->vec_func;
868646217ecSPeter Brune   B = snes->vec_rhs;
869fa9694d7SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
870d1adcc6fSPeter Brune   if (fas->downsmooth) {
871d1adcc6fSPeter Brune     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
872742fe5e2SPeter Brune     /* check convergence reason for the smoother */
873742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
874742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
875742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
876742fe5e2SPeter Brune       PetscFunctionReturn(0);
877742fe5e2SPeter Brune     }
878d28543b3SPeter Brune   } else if (snes->usegs && snes->ops->computegs) {
879794bee33SPeter Brune     if (fas->monitor) {
880794bee33SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
881794bee33SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
882d1adcc6fSPeter Brune       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
883eff52c0eSPeter Brune       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
884d1adcc6fSPeter Brune       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
885794bee33SPeter Brune     }
886646217ecSPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
887646217ecSPeter Brune       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
888cc05f883SPeter Brune       if (fas->monitor) {
889794bee33SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
890794bee33SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
891d1adcc6fSPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
892eff52c0eSPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
893d1adcc6fSPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
894794bee33SPeter Brune       }
895646217ecSPeter Brune     }
896c90fad12SPeter Brune   } else if (snes->pc) {
897c90fad12SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
8982f7ea302SPeter Brune     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
8992f7ea302SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
9002f7ea302SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
9012f7ea302SPeter Brune       PetscFunctionReturn(0);
9022f7ea302SPeter Brune     }
903fe6f9142SPeter Brune   }
904fa9694d7SPeter Brune   if (fas->next) {
905ee78dd50SPeter Brune     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
906794bee33SPeter Brune 
907c90fad12SPeter Brune     X_c  = fas->next->vec_sol;
908293a7e31SPeter Brune     Xo_c = fas->next->work[0];
909c90fad12SPeter Brune     F_c  = fas->next->vec_func;
910742fe5e2SPeter Brune     B_c  = fas->next->vec_rhs;
911efe1f98aSPeter Brune 
912efe1f98aSPeter Brune     /* inject the solution */
913efe1f98aSPeter Brune     if (fas->inject) {
914a3cb90a9SPeter Brune       ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
915efe1f98aSPeter Brune     } else {
916a3cb90a9SPeter Brune       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
917a3cb90a9SPeter Brune       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
918efe1f98aSPeter Brune     }
919293a7e31SPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
920293a7e31SPeter Brune 
921293a7e31SPeter Brune     /* restrict the defect */
922293a7e31SPeter Brune     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
923293a7e31SPeter Brune 
924c90fad12SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
925ee78dd50SPeter Brune     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
926c90fad12SPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
927742fe5e2SPeter Brune 
928293a7e31SPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
929c90fad12SPeter Brune 
930ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
931ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
932c90fad12SPeter Brune 
933c90fad12SPeter Brune     /* recurse to the next level */
934f5a6d4f9SBarry Smith     fas->next->vec_rhs = B_c;
935742fe5e2SPeter Brune     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
936742fe5e2SPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
937742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
938742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
939742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
940742fe5e2SPeter Brune       PetscFunctionReturn(0);
941742fe5e2SPeter Brune     }
942742fe5e2SPeter Brune     /* fas->next->vec_rhs = PETSC_NULL; */
943ee78dd50SPeter Brune 
944fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
945fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
946ee78dd50SPeter Brune     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X);CHKERRQ(ierr);
947293a7e31SPeter Brune   }
948742fe5e2SPeter Brune     /* up-smooth -- just update using the up-smoother */
949c90fad12SPeter Brune   if (fas->level != 0) {
950d1adcc6fSPeter Brune     if (fas->upsmooth) {
951d1adcc6fSPeter Brune       ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
952742fe5e2SPeter Brune       ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr);
953742fe5e2SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
954742fe5e2SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
955742fe5e2SPeter Brune         PetscFunctionReturn(0);
956742fe5e2SPeter Brune       }
957d28543b3SPeter Brune     } else if (snes->usegs && snes->ops->computegs) {
958794bee33SPeter Brune       if (fas->monitor) {
959794bee33SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
960794bee33SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
961d1adcc6fSPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
962eff52c0eSPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
963d1adcc6fSPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
964794bee33SPeter Brune       }
965646217ecSPeter Brune       for (k = 0; k < fas->max_down_it; k++) {
966646217ecSPeter Brune         ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
967794bee33SPeter Brune         if (fas->monitor) {
968794bee33SPeter Brune           ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
969794bee33SPeter Brune           ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
970d1adcc6fSPeter Brune           ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
971eff52c0eSPeter Brune           ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
972d1adcc6fSPeter Brune           ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
973794bee33SPeter Brune         }
974646217ecSPeter Brune       }
975c90fad12SPeter Brune     } else if (snes->pc) {
976c90fad12SPeter Brune       ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
9772f7ea302SPeter Brune       ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
9782f7ea302SPeter Brune       if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
9792f7ea302SPeter Brune         snes->reason = SNES_DIVERGED_INNER;
9802f7ea302SPeter Brune         PetscFunctionReturn(0);
9812f7ea302SPeter Brune       }
982c90fad12SPeter Brune     }
983fe6f9142SPeter Brune   }
984fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
985fa9694d7SPeter Brune 
986fa9694d7SPeter Brune   PetscFunctionReturn(0);
987421d9b32SPeter Brune }
988421d9b32SPeter Brune 
989421d9b32SPeter Brune #undef __FUNCT__
990421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
991421d9b32SPeter Brune 
992421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
993421d9b32SPeter Brune {
994fa9694d7SPeter Brune   PetscErrorCode ierr;
995fe6f9142SPeter Brune   PetscInt i, maxits;
996f5a6d4f9SBarry Smith   Vec X, B,F;
997fe6f9142SPeter Brune   PetscReal fnorm;
998421d9b32SPeter Brune   PetscFunctionBegin;
999fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
1000fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1001fa9694d7SPeter Brune   X = snes->vec_sol;
1002fe6f9142SPeter Brune   B = snes->vec_rhs;
1003f5a6d4f9SBarry Smith   F = snes->vec_func;
1004293a7e31SPeter Brune 
1005293a7e31SPeter Brune   /*norm setup */
1006fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1007fe6f9142SPeter Brune   snes->iter = 0;
1008fe6f9142SPeter Brune   snes->norm = 0.;
1009fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1010fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1011fe6f9142SPeter Brune   if (snes->domainerror) {
1012fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1013fe6f9142SPeter Brune     PetscFunctionReturn(0);
1014fe6f9142SPeter Brune   }
1015fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1016fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1017fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1018fe6f9142SPeter Brune   snes->norm = fnorm;
1019fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1020fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1021fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1022fe6f9142SPeter Brune 
1023fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
1024fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
1025fe6f9142SPeter Brune   /* test convergence */
1026fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1027fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
1028fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
1029fe6f9142SPeter Brune     /* Call general purpose update function */
1030646217ecSPeter Brune 
1031fe6f9142SPeter Brune     if (snes->ops->update) {
1032fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1033fe6f9142SPeter Brune     }
1034646217ecSPeter Brune     ierr = FASCycle_Private(snes, X);CHKERRQ(ierr);
1035742fe5e2SPeter Brune 
1036742fe5e2SPeter Brune     /* check for FAS cycle divergence */
1037742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
1038742fe5e2SPeter Brune       PetscFunctionReturn(0);
1039742fe5e2SPeter Brune     }
1040c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1041c90fad12SPeter Brune     /* Monitor convergence */
1042c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1043c90fad12SPeter Brune     snes->iter = i+1;
1044c90fad12SPeter Brune     snes->norm = fnorm;
1045c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1046c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
1047c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1048c90fad12SPeter Brune     /* Test for convergence */
1049c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1050c90fad12SPeter Brune     if (snes->reason) break;
1051fe6f9142SPeter Brune   }
1052fe6f9142SPeter Brune   if (i == maxits) {
1053fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1054fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1055fe6f9142SPeter Brune   }
1056421d9b32SPeter Brune   PetscFunctionReturn(0);
1057421d9b32SPeter Brune }
1058