xref: /petsc/src/snes/impls/fas/fas.c (revision 1fbfccc64754fe06f18ca1a892bc44eff051222f)
1421d9b32SPeter Brune /* Defines the basic SNES object */
26038b46bSPeter Brune #include <../src/snes/impls/fas/fasimpls.h>    /*I  "petscsnesfas.h"  I*/
3421d9b32SPeter Brune 
407144faaSPeter Brune const char *SNESFASTypes[] = {"MULTIPLICATIVE","ADDITIVE","SNESFASType","SNES_FAS",0};
507144faaSPeter Brune 
6421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
7421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
8421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
9421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
10421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
11421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
126273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *);
13421d9b32SPeter Brune 
14*1fbfccc6SJed Brown /*MC
15*1fbfccc6SJed Brown 
16*1fbfccc6SJed Brown SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
17*1fbfccc6SJed Brown 
18*1fbfccc6SJed Brown The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
19*1fbfccc6SJed Brown 
20*1fbfccc6SJed Brown Level: advanced
21*1fbfccc6SJed Brown 
22*1fbfccc6SJed Brown .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
23*1fbfccc6SJed Brown M*/
24421d9b32SPeter Brune 
25421d9b32SPeter Brune #undef __FUNCT__
26421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
27*1fbfccc6SJed Brown PETSC_EXTERN_C PetscErrorCode SNESCreate_FAS(SNES snes)
28421d9b32SPeter Brune {
29421d9b32SPeter Brune   SNES_FAS * fas;
30421d9b32SPeter Brune   PetscErrorCode ierr;
31421d9b32SPeter Brune 
32421d9b32SPeter Brune   PetscFunctionBegin;
33421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
34421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
35421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
36421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
37421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
38421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
39421d9b32SPeter Brune 
40ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
41ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
42ed020824SBarry Smith 
430e444f03SPeter Brune   snes->max_funcs = 30000;
440e444f03SPeter Brune   snes->max_its   = 10000;
450e444f03SPeter Brune 
46421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
47421d9b32SPeter Brune   snes->data                  = (void*) fas;
48421d9b32SPeter Brune   fas->level                  = 0;
49293a7e31SPeter Brune   fas->levels                 = 1;
50ee78dd50SPeter Brune   fas->n_cycles               = 1;
51ee78dd50SPeter Brune   fas->max_up_it              = 1;
52ee78dd50SPeter Brune   fas->max_down_it            = 1;
53ee78dd50SPeter Brune   fas->upsmooth               = PETSC_NULL;
54ee78dd50SPeter Brune   fas->downsmooth             = PETSC_NULL;
55421d9b32SPeter Brune   fas->next                   = PETSC_NULL;
566273346dSPeter Brune   fas->previous               = PETSC_NULL;
57421d9b32SPeter Brune   fas->interpolate            = PETSC_NULL;
58421d9b32SPeter Brune   fas->restrct                = PETSC_NULL;
59efe1f98aSPeter Brune   fas->inject                 = PETSC_NULL;
60cc05f883SPeter Brune   fas->monitor                = PETSC_NULL;
61cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
62ddebd997SPeter Brune   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
6322c1e704SPeter Brune   fas->linesearch_smooth      = PETSC_NULL;
64ddebd997SPeter Brune 
65421d9b32SPeter Brune   PetscFunctionReturn(0);
66421d9b32SPeter Brune }
67421d9b32SPeter Brune 
68421d9b32SPeter Brune #undef __FUNCT__
69421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
70c79ef259SPeter Brune /*@
712e8ce248SJed Brown    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
72c79ef259SPeter Brune 
73c79ef259SPeter Brune    Input Parameter:
742e8ce248SJed Brown .  snes - the nonlinear solver context
75c79ef259SPeter Brune 
76c79ef259SPeter Brune    Output parameter:
77c79ef259SPeter Brune .  levels - the number of levels
78c79ef259SPeter Brune 
79c79ef259SPeter Brune    Level: advanced
80c79ef259SPeter Brune 
81c79ef259SPeter Brune .keywords: MG, get, levels, multigrid
82c79ef259SPeter Brune 
83c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels()
84c79ef259SPeter Brune @*/
85421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
86421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
87421d9b32SPeter Brune   PetscFunctionBegin;
88ee1fd11aSPeter Brune   *levels = fas->levels;
89421d9b32SPeter Brune   PetscFunctionReturn(0);
90421d9b32SPeter Brune }
91421d9b32SPeter Brune 
92421d9b32SPeter Brune #undef __FUNCT__
93646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles"
94c79ef259SPeter Brune /*@
952e8ce248SJed Brown    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
96c79ef259SPeter Brune    complicated cycling.
97c79ef259SPeter Brune 
98c79ef259SPeter Brune    Logically Collective on SNES
99c79ef259SPeter Brune 
100c79ef259SPeter Brune    Input Parameters:
101c79ef259SPeter Brune +  snes   - the multigrid context
102c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
103c79ef259SPeter Brune 
104c79ef259SPeter Brune    Options Database Key:
105c79ef259SPeter Brune $  -snes_fas_cycles 1 or 2
106c79ef259SPeter Brune 
107c79ef259SPeter Brune    Level: advanced
108c79ef259SPeter Brune 
109c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
110c79ef259SPeter Brune 
111c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel()
112c79ef259SPeter Brune @*/
113646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
114646217ecSPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
115646217ecSPeter Brune   PetscErrorCode ierr;
116646217ecSPeter Brune   PetscFunctionBegin;
117646217ecSPeter Brune   fas->n_cycles = cycles;
118646217ecSPeter Brune   if (fas->next) {
119646217ecSPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
120646217ecSPeter Brune   }
121646217ecSPeter Brune   PetscFunctionReturn(0);
122646217ecSPeter Brune }
123646217ecSPeter Brune 
124eff52c0eSPeter Brune #undef __FUNCT__
125c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel"
126c79ef259SPeter Brune /*@
127722262beSPeter Brune    SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level.
128c79ef259SPeter Brune 
129c79ef259SPeter Brune    Logically Collective on SNES
130c79ef259SPeter Brune 
131c79ef259SPeter Brune    Input Parameters:
132c79ef259SPeter Brune +  snes   - the multigrid context
133c79ef259SPeter Brune .  level  - the level to set the number of cycles on
134c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
135c79ef259SPeter Brune 
136c79ef259SPeter Brune    Level: advanced
137c79ef259SPeter Brune 
138c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
139c79ef259SPeter Brune 
140c79ef259SPeter Brune .seealso: SNESFASSetCycles()
141c79ef259SPeter Brune @*/
142c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) {
143c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
144c79ef259SPeter Brune   PetscInt top_level = fas->level,i;
145c79ef259SPeter Brune 
146c79ef259SPeter Brune   PetscFunctionBegin;
1472e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D cycle count from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
148c79ef259SPeter Brune   /* get to the correct level */
149c79ef259SPeter Brune   for (i = fas->level; i > level; i--) {
150c79ef259SPeter Brune     fas = (SNES_FAS *)fas->next->data;
151c79ef259SPeter Brune   }
1522e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
153c79ef259SPeter Brune   fas->n_cycles = cycles;
154c79ef259SPeter Brune   PetscFunctionReturn(0);
155c79ef259SPeter Brune }
156c79ef259SPeter Brune 
157c79ef259SPeter Brune #undef __FUNCT__
158eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS"
159aeed3662SMatthew G Knepley /*@C
160c79ef259SPeter Brune    SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used.
161c79ef259SPeter Brune    Use SNESFASSetGSOnLevel() for more complicated staging of smoothers
162c79ef259SPeter Brune    and nonlinear preconditioners.
163c79ef259SPeter Brune 
164c79ef259SPeter Brune    Logically Collective on SNES
165c79ef259SPeter Brune 
166c79ef259SPeter Brune    Input Parameters:
167c79ef259SPeter Brune +  snes    - the multigrid context
168c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
169c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
170c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
171c79ef259SPeter Brune 
172c79ef259SPeter Brune    Level: advanced
173c79ef259SPeter Brune 
174c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid
175c79ef259SPeter Brune 
176c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel()
177c79ef259SPeter Brune @*/
178eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
179eff52c0eSPeter Brune   PetscErrorCode ierr = 0;
180eff52c0eSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
181eff52c0eSPeter Brune   PetscFunctionBegin;
182eff52c0eSPeter Brune 
183eff52c0eSPeter Brune   if (gsfunc) {
184eff52c0eSPeter Brune     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
185eff52c0eSPeter Brune     /* push the provided GS up the tree */
186eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
1876cab3a1bSJed Brown   } else {
1886cab3a1bSJed Brown     ierr = SNESGetGS(snes,&gsfunc,&ctx);CHKERRQ(ierr);
1896cab3a1bSJed Brown     if (gsfunc) {
190eff52c0eSPeter Brune       /* assume that the user has set the GS solver at this level */
191eff52c0eSPeter Brune       if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
192eff52c0eSPeter Brune     } else if (use_gs) {
1932e8ce248SJed Brown       SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %D", fas->level);
194eff52c0eSPeter Brune     }
1956cab3a1bSJed Brown   }
196eff52c0eSPeter Brune   PetscFunctionReturn(0);
197eff52c0eSPeter Brune }
198eff52c0eSPeter Brune 
199eff52c0eSPeter Brune #undef __FUNCT__
200eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel"
201aeed3662SMatthew G Knepley /*@C
202c79ef259SPeter Brune    SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level.
203c79ef259SPeter Brune 
204c79ef259SPeter Brune    Logically Collective on SNES
205c79ef259SPeter Brune 
206c79ef259SPeter Brune    Input Parameters:
207c79ef259SPeter Brune +  snes    - the multigrid context
208c79ef259SPeter Brune .  level   - the level to set the nonlinear smoother on
209c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
210c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
211c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
212c79ef259SPeter Brune 
213c79ef259SPeter Brune    Level: advanced
214c79ef259SPeter Brune 
215c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
216c79ef259SPeter Brune 
217c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS()
218c79ef259SPeter Brune @*/
219eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
220eff52c0eSPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
221eff52c0eSPeter Brune   PetscErrorCode ierr;
222eff52c0eSPeter Brune   PetscInt       top_level = fas->level,i;
223eff52c0eSPeter Brune   SNES           cur_snes = snes;
224eff52c0eSPeter Brune   PetscFunctionBegin;
2252e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D GS from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
226eff52c0eSPeter Brune   /* get to the correct level */
227eff52c0eSPeter Brune   for (i = fas->level; i > level; i--) {
228eff52c0eSPeter Brune     fas = (SNES_FAS *)fas->next->data;
229eff52c0eSPeter Brune     cur_snes = fas->next;
230eff52c0eSPeter Brune   }
2312e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
232eff52c0eSPeter Brune   if (gsfunc) {
2336273346dSPeter Brune     ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr);
234eff52c0eSPeter Brune   }
235eff52c0eSPeter Brune   PetscFunctionReturn(0);
236eff52c0eSPeter Brune }
237eff52c0eSPeter Brune 
238646217ecSPeter Brune #undef __FUNCT__
239421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
240c79ef259SPeter Brune /*@
241c79ef259SPeter Brune    SNESFASGetSNES - Gets the SNES corresponding to a particular
242c79ef259SPeter Brune    level of the FAS hierarchy.
243c79ef259SPeter Brune 
244c79ef259SPeter Brune    Input Parameters:
245c79ef259SPeter Brune +  snes    - the multigrid context
246c79ef259SPeter Brune    level   - the level to get
247c79ef259SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
248c79ef259SPeter Brune 
249c79ef259SPeter Brune    Level: advanced
250c79ef259SPeter Brune 
251c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
252c79ef259SPeter Brune 
253c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
254c79ef259SPeter Brune @*/
255421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes,PetscInt level,SNES *lsnes) {
256421d9b32SPeter Brune   SNES_FAS *fas = (SNES_FAS*)snes->data;
257421d9b32SPeter Brune   PetscInt i;
2582e8ce248SJed Brown 
259421d9b32SPeter Brune   PetscFunctionBegin;
2602e8ce248SJed Brown   if (level > fas->levels-1) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
2612e8ce248SJed Brown   if (fas->level < level) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS on level %D, must call from SNES at least as fine as level",level,fas->level);
262421d9b32SPeter Brune   *lsnes = snes;
263421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
264421d9b32SPeter Brune     *lsnes = fas->next;
265421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
266421d9b32SPeter Brune   }
2672e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
268421d9b32SPeter Brune   PetscFunctionReturn(0);
269421d9b32SPeter Brune }
270421d9b32SPeter Brune 
271421d9b32SPeter Brune #undef __FUNCT__
27207144faaSPeter Brune #define __FUNCT__ "SNESFASSetType"
27307144faaSPeter Brune /*@
27407144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS.
27507144faaSPeter Brune 
276*1fbfccc6SJed Brown Logically Collective
27707144faaSPeter Brune 
278*1fbfccc6SJed Brown Input Arguments:
279*1fbfccc6SJed Brown snes - nonlinear solver
280*1fbfccc6SJed Brown fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE
281*1fbfccc6SJed Brown 
282*1fbfccc6SJed Brown .seealso: PCMGSetType()
28307144faaSPeter Brune @*/
28407144faaSPeter Brune PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
28507144faaSPeter Brune {
28607144faaSPeter Brune   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
28707144faaSPeter Brune 
28807144faaSPeter Brune   PetscFunctionBegin;
28907144faaSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
29007144faaSPeter Brune   PetscValidLogicalCollectiveEnum(snes,fastype,2);
29107144faaSPeter Brune   fas->fastype = fastype;
29207144faaSPeter Brune   PetscFunctionReturn(0);
29307144faaSPeter Brune }
29407144faaSPeter Brune 
29507144faaSPeter Brune #undef __FUNCT__
296421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
297aeed3662SMatthew G Knepley /*@C
298c79ef259SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
299c79ef259SPeter Brune    Must be called before any other FAS routine.
300c79ef259SPeter Brune 
301c79ef259SPeter Brune    Input Parameters:
302c79ef259SPeter Brune +  snes   - the snes context
303c79ef259SPeter Brune .  levels - the number of levels
304c79ef259SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
305c79ef259SPeter Brune             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
306c79ef259SPeter Brune             Fortran.
307c79ef259SPeter Brune 
308c79ef259SPeter Brune    Level: intermediate
309c79ef259SPeter Brune 
310c79ef259SPeter Brune    Notes:
311c79ef259SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
312c79ef259SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
313c79ef259SPeter Brune 
314c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid
315c79ef259SPeter Brune 
316c79ef259SPeter Brune .seealso: SNESFASGetLevels()
317c79ef259SPeter Brune @*/
318421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
319421d9b32SPeter Brune   PetscErrorCode ierr;
320421d9b32SPeter Brune   PetscInt i;
321421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
3226273346dSPeter Brune   SNES prevsnes;
323421d9b32SPeter Brune   MPI_Comm comm;
324421d9b32SPeter Brune   PetscFunctionBegin;
325ee1fd11aSPeter Brune   comm = ((PetscObject)snes)->comm;
326ee1fd11aSPeter Brune   if (levels == fas->levels) {
327ee1fd11aSPeter Brune     if (!comms)
328ee1fd11aSPeter Brune       PetscFunctionReturn(0);
329ee1fd11aSPeter Brune   }
330421d9b32SPeter Brune   /* user has changed the number of levels; reset */
331421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
332421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
333421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
334ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
3356273346dSPeter Brune   fas->previous = PETSC_NULL;
3366273346dSPeter Brune   prevsnes = snes;
337421d9b32SPeter Brune   /* setup the finest level */
338421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
339421d9b32SPeter Brune     if (comms) comm = comms[i];
340421d9b32SPeter Brune     fas->level = i;
341421d9b32SPeter Brune     fas->levels = levels;
342ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
343421d9b32SPeter Brune     if (i > 0) {
344421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
3454a6b3fb8SBarry Smith       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
346421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
347794bee33SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
3486273346dSPeter Brune       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
3496273346dSPeter Brune       prevsnes = fas->next;
3506273346dSPeter Brune       fas = (SNES_FAS *)prevsnes->data;
351421d9b32SPeter Brune     }
352421d9b32SPeter Brune   }
353421d9b32SPeter Brune   PetscFunctionReturn(0);
354421d9b32SPeter Brune }
355421d9b32SPeter Brune 
356421d9b32SPeter Brune #undef __FUNCT__
357c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp"
358c79ef259SPeter Brune /*@
359c79ef259SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
360c79ef259SPeter Brune    use on all levels.
361c79ef259SPeter Brune 
362fde0ff24SPeter Brune    Logically Collective on SNES
363c79ef259SPeter Brune 
364c79ef259SPeter Brune    Input Parameters:
365c79ef259SPeter Brune +  snes - the multigrid context
366c79ef259SPeter Brune -  n    - the number of smoothing steps
367c79ef259SPeter Brune 
368c79ef259SPeter Brune    Options Database Key:
369d28543b3SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
370c79ef259SPeter Brune 
371c79ef259SPeter Brune    Level: advanced
372c79ef259SPeter Brune 
373c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
374c79ef259SPeter Brune 
375c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
376c79ef259SPeter Brune @*/
377c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
378c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
379c79ef259SPeter Brune   PetscErrorCode ierr = 0;
380c79ef259SPeter Brune   PetscFunctionBegin;
381c79ef259SPeter Brune   fas->max_up_it = n;
382c79ef259SPeter Brune   if (fas->next) {
383c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
384c79ef259SPeter Brune   }
385c79ef259SPeter Brune   PetscFunctionReturn(0);
386c79ef259SPeter Brune }
387c79ef259SPeter Brune 
388c79ef259SPeter Brune #undef __FUNCT__
389c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown"
390c79ef259SPeter Brune /*@
391c79ef259SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
392c79ef259SPeter Brune    use on all levels.
393c79ef259SPeter Brune 
394fde0ff24SPeter Brune    Logically Collective on SNES
395c79ef259SPeter Brune 
396c79ef259SPeter Brune    Input Parameters:
397c79ef259SPeter Brune +  snes - the multigrid context
398c79ef259SPeter Brune -  n    - the number of smoothing steps
399c79ef259SPeter Brune 
400c79ef259SPeter Brune    Options Database Key:
401d28543b3SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
402c79ef259SPeter Brune 
403c79ef259SPeter Brune    Level: advanced
404c79ef259SPeter Brune 
405c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
406c79ef259SPeter Brune 
407c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
408c79ef259SPeter Brune @*/
409c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
410c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
411c79ef259SPeter Brune   PetscErrorCode ierr = 0;
412c79ef259SPeter Brune   PetscFunctionBegin;
413c79ef259SPeter Brune   fas->max_down_it = n;
414c79ef259SPeter Brune   if (fas->next) {
415c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
416c79ef259SPeter Brune   }
417c79ef259SPeter Brune   PetscFunctionReturn(0);
418c79ef259SPeter Brune }
419c79ef259SPeter Brune 
420c79ef259SPeter Brune #undef __FUNCT__
421421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
422c79ef259SPeter Brune /*@
423c79ef259SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
424c79ef259SPeter Brune    interpolation from l-1 to the lth level
425c79ef259SPeter Brune 
426c79ef259SPeter Brune    Input Parameters:
427c79ef259SPeter Brune +  snes      - the multigrid context
428c79ef259SPeter Brune .  mat       - the interpolation operator
429c79ef259SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
430c79ef259SPeter Brune 
431c79ef259SPeter Brune    Level: advanced
432c79ef259SPeter Brune 
433c79ef259SPeter Brune    Notes:
434c79ef259SPeter Brune           Usually this is the same matrix used also to set the restriction
435c79ef259SPeter Brune     for the same level.
436c79ef259SPeter Brune 
437c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
438c79ef259SPeter Brune     out from the matrix size which one it is.
439c79ef259SPeter Brune 
440c79ef259SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
441c79ef259SPeter Brune 
442bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
443c79ef259SPeter Brune @*/
444421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
445421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
446d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
447bccf9bb3SJed Brown   PetscErrorCode ierr;
448421d9b32SPeter Brune 
449421d9b32SPeter Brune   PetscFunctionBegin;
4502e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D interpolation from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
451421d9b32SPeter Brune   /* get to the correct level */
452d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
453421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
454421d9b32SPeter Brune   }
4552e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
456bccf9bb3SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
457bd4e12b0SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
458421d9b32SPeter Brune   fas->interpolate = mat;
459421d9b32SPeter Brune   PetscFunctionReturn(0);
460421d9b32SPeter Brune }
461421d9b32SPeter Brune 
462421d9b32SPeter Brune #undef __FUNCT__
463421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
464c79ef259SPeter Brune /*@
465c79ef259SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
466c79ef259SPeter Brune    from level l to l-1.
467c79ef259SPeter Brune 
468c79ef259SPeter Brune    Input Parameters:
469c79ef259SPeter Brune +  snes  - the multigrid context
470c79ef259SPeter Brune .  mat   - the restriction matrix
471c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
472c79ef259SPeter Brune 
473c79ef259SPeter Brune    Level: advanced
474c79ef259SPeter Brune 
475c79ef259SPeter Brune    Notes:
476c79ef259SPeter Brune           Usually this is the same matrix used also to set the interpolation
477c79ef259SPeter Brune     for the same level.
478c79ef259SPeter Brune 
479c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
480c79ef259SPeter Brune     out from the matrix size which one it is.
481c79ef259SPeter Brune 
482fde0ff24SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
483c79ef259SPeter Brune     is used.
484c79ef259SPeter Brune 
485c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
486c79ef259SPeter Brune 
487c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
488c79ef259SPeter Brune @*/
489421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
490421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
491d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
492bccf9bb3SJed Brown   PetscErrorCode ierr;
493421d9b32SPeter Brune 
494421d9b32SPeter Brune   PetscFunctionBegin;
4952e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D restriction from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
496421d9b32SPeter Brune   /* get to the correct level */
497d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
498421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
499421d9b32SPeter Brune   }
5002e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
501bccf9bb3SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
502bd4e12b0SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
503421d9b32SPeter Brune   fas->restrct = mat;
504421d9b32SPeter Brune   PetscFunctionReturn(0);
505421d9b32SPeter Brune }
506421d9b32SPeter Brune 
507421d9b32SPeter Brune #undef __FUNCT__
508421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
509c79ef259SPeter Brune /*@
510c79ef259SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
511c79ef259SPeter Brune    operator from level l to l-1.
512c79ef259SPeter Brune 
513c79ef259SPeter Brune    Input Parameters:
514c79ef259SPeter Brune +  snes   - the multigrid context
515c79ef259SPeter Brune .  rscale - the restriction scaling
516c79ef259SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
517c79ef259SPeter Brune 
518c79ef259SPeter Brune    Level: advanced
519c79ef259SPeter Brune 
520c79ef259SPeter Brune    Notes:
521c79ef259SPeter Brune          This is only used in the case that the injection is not set.
522c79ef259SPeter Brune 
523c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
524c79ef259SPeter Brune 
525c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
526c79ef259SPeter Brune @*/
527421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
528421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
529d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
530bccf9bb3SJed Brown   PetscErrorCode ierr;
531421d9b32SPeter Brune 
532421d9b32SPeter Brune   PetscFunctionBegin;
5332e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D rscale from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
534421d9b32SPeter Brune   /* get to the correct level */
535d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
536421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
537421d9b32SPeter Brune   }
5382e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
539bccf9bb3SJed Brown   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
540bd4e12b0SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
541421d9b32SPeter Brune   fas->rscale = rscale;
542421d9b32SPeter Brune   PetscFunctionReturn(0);
543421d9b32SPeter Brune }
544421d9b32SPeter Brune 
545efe1f98aSPeter Brune 
546efe1f98aSPeter Brune #undef __FUNCT__
547efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection"
548c79ef259SPeter Brune /*@
549c79ef259SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
550c79ef259SPeter Brune    from level l to l-1.
551c79ef259SPeter Brune 
552c79ef259SPeter Brune    Input Parameters:
553c79ef259SPeter Brune +  snes  - the multigrid context
554c79ef259SPeter Brune .  mat   - the restriction matrix
555c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
556c79ef259SPeter Brune 
557c79ef259SPeter Brune    Level: advanced
558c79ef259SPeter Brune 
559c79ef259SPeter Brune    Notes:
560c79ef259SPeter Brune          If you do not set this, the restriction and rscale is used to
561c79ef259SPeter Brune    project the solution instead.
562c79ef259SPeter Brune 
563c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
564c79ef259SPeter Brune 
565c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
566c79ef259SPeter Brune @*/
567efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
568efe1f98aSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
569efe1f98aSPeter Brune   PetscInt top_level = fas->level,i;
570bccf9bb3SJed Brown   PetscErrorCode ierr;
571efe1f98aSPeter Brune 
572efe1f98aSPeter Brune   PetscFunctionBegin;
5732e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D injection from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
574efe1f98aSPeter Brune   /* get to the correct level */
575efe1f98aSPeter Brune   for (i = fas->level; i > level; i--) {
576efe1f98aSPeter Brune     fas = (SNES_FAS *)fas->next->data;
577efe1f98aSPeter Brune   }
5782e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
579bccf9bb3SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
580bd4e12b0SJed Brown   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
581efe1f98aSPeter Brune   fas->inject = mat;
582efe1f98aSPeter Brune   PetscFunctionReturn(0);
583efe1f98aSPeter Brune }
584efe1f98aSPeter Brune 
585421d9b32SPeter Brune #undef __FUNCT__
586421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
587421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
588421d9b32SPeter Brune {
58977df8cc4SPeter Brune   PetscErrorCode ierr = 0;
590421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
591421d9b32SPeter Brune 
592421d9b32SPeter Brune   PetscFunctionBegin;
593bccf9bb3SJed Brown   ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
594bccf9bb3SJed Brown   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
5953dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
596bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
597bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
598bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
5993dccd265SPeter Brune 
6003dccd265SPeter Brune   if (fas->upsmooth)   ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
6013dccd265SPeter Brune   if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
602742fe5e2SPeter Brune   if (fas->next)       ierr = SNESReset(fas->next);CHKERRQ(ierr);
60322c1e704SPeter Brune 
6046188f407SPeter Brune   ierr = PetscLineSearchDestroy(&fas->linesearch_smooth);CHKERRQ(ierr);
60522c1e704SPeter Brune 
606421d9b32SPeter Brune   PetscFunctionReturn(0);
607421d9b32SPeter Brune }
608421d9b32SPeter Brune 
609421d9b32SPeter Brune #undef __FUNCT__
610421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
611421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
612421d9b32SPeter Brune {
613421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
614742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
615421d9b32SPeter Brune 
616421d9b32SPeter Brune   PetscFunctionBegin;
617421d9b32SPeter Brune   /* recursively resets and then destroys */
61879d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
619421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
620421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
621ee1fd11aSPeter Brune 
622421d9b32SPeter Brune   PetscFunctionReturn(0);
623421d9b32SPeter Brune }
624421d9b32SPeter Brune 
625421d9b32SPeter Brune #undef __FUNCT__
626421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
627421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
628421d9b32SPeter Brune {
62948bfdf8aSPeter Brune   SNES_FAS                *fas = (SNES_FAS *) snes->data;
630421d9b32SPeter Brune   PetscErrorCode          ierr;
63122c1e704SPeter Brune   const char              *optionsprefix;
632efe1f98aSPeter Brune   VecScatter              injscatter;
633d1adcc6fSPeter Brune   PetscInt                dm_levels;
6343dccd265SPeter Brune   Vec                     vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
635d1adcc6fSPeter Brune 
636421d9b32SPeter Brune   PetscFunctionBegin;
637eff52c0eSPeter Brune 
638cc05f883SPeter Brune   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
639d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
640d1adcc6fSPeter Brune     dm_levels++;
641cc05f883SPeter Brune     if (dm_levels > fas->levels) {
6423dccd265SPeter Brune 
6432e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
6443dccd265SPeter Brune       vec_sol = snes->vec_sol;
6453dccd265SPeter Brune       vec_func = snes->vec_func;
6463dccd265SPeter Brune       vec_sol_update = snes->vec_sol_update;
6473dccd265SPeter Brune       vec_rhs = snes->vec_rhs;
6483dccd265SPeter Brune       snes->vec_sol = PETSC_NULL;
6493dccd265SPeter Brune       snes->vec_func = PETSC_NULL;
6503dccd265SPeter Brune       snes->vec_sol_update = PETSC_NULL;
6513dccd265SPeter Brune       snes->vec_rhs = PETSC_NULL;
6523dccd265SPeter Brune 
6533dccd265SPeter Brune       /* reset the number of levels */
654d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
655cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
6563dccd265SPeter Brune 
6573dccd265SPeter Brune       snes->vec_sol = vec_sol;
6583dccd265SPeter Brune       snes->vec_func = vec_func;
6593dccd265SPeter Brune       snes->vec_rhs = vec_rhs;
6603dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
661d1adcc6fSPeter Brune     }
662d1adcc6fSPeter Brune   }
663d1adcc6fSPeter Brune 
6643dccd265SPeter Brune   if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
6653dccd265SPeter Brune 
66607144faaSPeter Brune   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
667fde0ff24SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
66807144faaSPeter Brune   } else {
669ddebd997SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
67007144faaSPeter Brune   }
671cc05f883SPeter Brune 
67279d9a41aSPeter Brune   if (snes->dm) {
67379d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
67479d9a41aSPeter Brune     if (fas->next) {
67579d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
67679d9a41aSPeter Brune       if (!fas->next->dm) {
67779d9a41aSPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
67879d9a41aSPeter Brune         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
67979d9a41aSPeter Brune       }
68079d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
68179d9a41aSPeter Brune       if (!fas->interpolate) {
68279d9a41aSPeter Brune         ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
683bccf9bb3SJed Brown         if (!fas->restrct) {
684bccf9bb3SJed Brown           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
68579d9a41aSPeter Brune           fas->restrct = fas->interpolate;
68679d9a41aSPeter Brune         }
687bccf9bb3SJed Brown       }
68879d9a41aSPeter Brune       /* set the injection from the DM */
68979d9a41aSPeter Brune       if (!fas->inject) {
69079d9a41aSPeter Brune         ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
69179d9a41aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
69279d9a41aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
69379d9a41aSPeter Brune       }
69479d9a41aSPeter Brune     }
69579d9a41aSPeter Brune     /* set the DMs of the pre and post-smoothers here */
69679d9a41aSPeter Brune     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
69779d9a41aSPeter Brune     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
69879d9a41aSPeter Brune   }
69979d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
70079d9a41aSPeter Brune 
70179d9a41aSPeter Brune  if (fas->next) {
70279d9a41aSPeter Brune     if (fas->galerkin) {
70379d9a41aSPeter Brune       ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr);
70479d9a41aSPeter Brune     }
70579d9a41aSPeter Brune   }
70679d9a41aSPeter Brune 
70779d9a41aSPeter Brune   if (fas->next) {
70879d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
709938e4a01SJed Brown     if (!fas->next->vec_sol) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_sol);CHKERRQ(ierr);}
710938e4a01SJed Brown     if (!fas->next->vec_rhs) {ierr = SNESFASCreateCoarseVec(snes,&fas->next->vec_rhs);CHKERRQ(ierr);}
71179d9a41aSPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
71279d9a41aSPeter Brune   }
71379d9a41aSPeter Brune 
7146cab3a1bSJed Brown   /* setup the pre and post smoothers */
7156cab3a1bSJed Brown   if (fas->upsmooth) {ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);}
7166cab3a1bSJed Brown   if (fas->downsmooth) {ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);}
717cc05f883SPeter Brune 
71822c1e704SPeter Brune   /* if the pre and post smoothers don't exist, set up line searches in their place */
71922c1e704SPeter Brune   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
72022c1e704SPeter Brune   if (!fas->upsmooth || !fas->downsmooth) {
7216188f407SPeter Brune     ierr = PetscLineSearchCreate(((PetscObject)snes)->comm, &fas->linesearch_smooth);CHKERRQ(ierr);
7226188f407SPeter Brune     ierr = PetscLineSearchSetSNES(fas->linesearch_smooth, snes);CHKERRQ(ierr);
7236188f407SPeter Brune     ierr = PetscLineSearchSetType(fas->linesearch_smooth, PETSCLINESEARCHL2);CHKERRQ(ierr);
7246188f407SPeter Brune     ierr = PetscLineSearchAppendOptionsPrefix(fas->linesearch_smooth, "fas_");CHKERRQ(ierr);
7256188f407SPeter Brune     ierr = PetscLineSearchAppendOptionsPrefix(fas->linesearch_smooth, optionsprefix);CHKERRQ(ierr);
7266188f407SPeter Brune     ierr = PetscLineSearchSetFromOptions(fas->linesearch_smooth);CHKERRQ(ierr);
72722c1e704SPeter Brune   }
72822c1e704SPeter Brune 
7296273346dSPeter Brune   /* setup FAS work vectors */
7306273346dSPeter Brune   if (fas->galerkin) {
7316273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
7326273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
7336273346dSPeter Brune   }
734421d9b32SPeter Brune   PetscFunctionReturn(0);
735421d9b32SPeter Brune }
736421d9b32SPeter Brune 
737421d9b32SPeter Brune #undef __FUNCT__
738421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
739421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
740421d9b32SPeter Brune {
741ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
742ee78dd50SPeter Brune   PetscInt       levels = 1;
743fde0ff24SPeter Brune   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg;
744421d9b32SPeter Brune   PetscErrorCode ierr;
745ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
74607144faaSPeter Brune   SNESFASType    fastype;
747fde0ff24SPeter Brune   const char     *optionsprefix;
748fde0ff24SPeter Brune   const char     *prefix;
7499e764e56SPeter Brune   PetscLineSearch linesearch;
750421d9b32SPeter Brune 
751421d9b32SPeter Brune   PetscFunctionBegin;
752c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
753ee78dd50SPeter Brune 
754ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
755ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
756ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
757c732cbdbSBarry Smith     if (!flg && snes->dm) {
758c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
759c732cbdbSBarry Smith       levels++;
760d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
761c732cbdbSBarry Smith     }
762ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
763ee78dd50SPeter Brune   }
76407144faaSPeter Brune   fastype = fas->fastype;
76507144faaSPeter Brune   ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
76607144faaSPeter Brune   if (flg) {
76707144faaSPeter Brune     ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
76807144faaSPeter Brune   }
769ee78dd50SPeter Brune 
770fde0ff24SPeter Brune   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
771fde0ff24SPeter Brune 
772fde0ff24SPeter Brune   /* smoother setup options */
773fde0ff24SPeter Brune   ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr);
774fde0ff24SPeter Brune   ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr);
775fde0ff24SPeter Brune   ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr);
776fde0ff24SPeter Brune   if (fas->level == 0) {
777fde0ff24SPeter Brune     ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr);
778fde0ff24SPeter Brune   }
779ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
780fde0ff24SPeter Brune 
781d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
782ee78dd50SPeter Brune 
7836273346dSPeter Brune   ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr);
7846273346dSPeter Brune 
785646217ecSPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
786ee78dd50SPeter Brune 
787ee78dd50SPeter Brune 
788421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
7898cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
790eff52c0eSPeter Brune 
791fde0ff24SPeter Brune   if ((!fas->downsmooth) && smoothcoarseflg) {
792fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
793fde0ff24SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
794fde0ff24SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
795fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
796fde0ff24SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
797fde0ff24SPeter Brune   }
798fde0ff24SPeter Brune 
799fde0ff24SPeter Brune   if ((!fas->downsmooth) && smoothdownflg) {
8008cc86e31SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
8018cc86e31SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
8028cc86e31SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
803fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr);
8048cc86e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
8058cc86e31SPeter Brune   }
8068cc86e31SPeter Brune 
807fde0ff24SPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) {
80867339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
809ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
81067339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
811fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr);
812293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
813ee78dd50SPeter Brune   }
814fde0ff24SPeter Brune 
815fde0ff24SPeter Brune   if ((!fas->downsmooth) && smoothflg) {
816fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
817fde0ff24SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
818fde0ff24SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
819fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr);
820fde0ff24SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
821fde0ff24SPeter Brune   }
822fde0ff24SPeter Brune 
823fde0ff24SPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) {
824fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
825fde0ff24SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
826fde0ff24SPeter Brune     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
827fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr);
828fde0ff24SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
829fde0ff24SPeter Brune   }
830fde0ff24SPeter Brune 
8311a266240SBarry Smith   if (fas->upsmooth) {
832fde0ff24SPeter Brune     ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr);
8331a266240SBarry Smith   }
8341a266240SBarry Smith 
8351a266240SBarry Smith   if (fas->downsmooth) {
836fde0ff24SPeter Brune     ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr);
8371a266240SBarry Smith   }
838ee78dd50SPeter Brune 
839742fe5e2SPeter Brune   if (fas->level != fas->levels - 1) {
840fde0ff24SPeter Brune     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr);
841742fe5e2SPeter Brune   }
842742fe5e2SPeter Brune 
843ce11c761SPeter Brune   /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */
844ce11c761SPeter Brune   if (!fas->downsmooth) {
84593dfaa4eSPeter Brune     ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
8465d115551SPeter Brune     ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
847ce11c761SPeter Brune     if (fas->level == 0) {
84879d9a41aSPeter Brune       ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
849ce11c761SPeter Brune     }
850ce11c761SPeter Brune   }
851ce11c761SPeter Brune 
852ce11c761SPeter Brune   if (!fas->upsmooth) {
85393dfaa4eSPeter Brune     ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
8545d115551SPeter Brune     ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
855ce11c761SPeter Brune   }
856ce11c761SPeter Brune 
857ee78dd50SPeter Brune   if (monflg) {
858646217ecSPeter Brune     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
859794bee33SPeter Brune     /* set the monitors for the upsmoother and downsmoother */
8602f7ea302SPeter Brune     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
861742fe5e2SPeter Brune     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
862293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
863293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
864d28543b3SPeter Brune   } else {
865d28543b3SPeter Brune     /* unset the monitors on the coarse levels */
866d28543b3SPeter Brune     if (fas->level != fas->levels - 1) {
867d28543b3SPeter Brune       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
868d28543b3SPeter Brune     }
869ee78dd50SPeter Brune   }
870ee78dd50SPeter Brune 
8719e764e56SPeter Brune 
8729e764e56SPeter Brune   /* set up the default line search for coarse grid corrections */
8739e764e56SPeter Brune   if (fas->fastype == SNES_FAS_ADDITIVE) {
8749e764e56SPeter Brune     if (!snes->linesearch) {
8759e764e56SPeter Brune       ierr = SNESGetPetscLineSearch(snes, &linesearch);CHKERRQ(ierr);
8769e764e56SPeter Brune       ierr = PetscLineSearchSetType(linesearch, PETSCLINESEARCHL2);CHKERRQ(ierr);
8779e764e56SPeter Brune     }
8789e764e56SPeter Brune   }
8799e764e56SPeter Brune 
880ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
881d28543b3SPeter Brune   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
882421d9b32SPeter Brune   PetscFunctionReturn(0);
883421d9b32SPeter Brune }
884421d9b32SPeter Brune 
885421d9b32SPeter Brune #undef __FUNCT__
886421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
887421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
888421d9b32SPeter Brune {
889421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
890421d9b32SPeter Brune   PetscBool      iascii;
891421d9b32SPeter Brune   PetscErrorCode ierr;
892421d9b32SPeter Brune 
893421d9b32SPeter Brune   PetscFunctionBegin;
894421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
895421d9b32SPeter Brune   if (iascii) {
8962e8ce248SJed Brown     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n",  fas->levels);CHKERRQ(ierr);
897421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
898bd4e12b0SJed Brown     ierr = PetscViewerASCIIPrintf(viewer, "level: %D\n",  fas->level);CHKERRQ(ierr);
899ee78dd50SPeter Brune     if (fas->upsmooth) {
90039a4b097SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
901421d9b32SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);
902ee78dd50SPeter Brune       ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
903421d9b32SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);
904421d9b32SPeter Brune     } else {
905ee78dd50SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
906421d9b32SPeter Brune     }
907ee78dd50SPeter Brune     if (fas->downsmooth) {
90839a4b097SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
909421d9b32SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);
910ee78dd50SPeter Brune       ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
911421d9b32SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);
912421d9b32SPeter Brune     } else {
913ee78dd50SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
914421d9b32SPeter Brune     }
915421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
916421d9b32SPeter Brune   } else {
917421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
918421d9b32SPeter Brune   }
919421d9b32SPeter Brune   PetscFunctionReturn(0);
920421d9b32SPeter Brune }
921421d9b32SPeter Brune 
922421d9b32SPeter Brune #undef __FUNCT__
92339bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth"
92439bd7f45SPeter Brune /*
92539bd7f45SPeter Brune Defines the action of the downsmoother
92639bd7f45SPeter Brune  */
92739bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
92839bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
92922c1e704SPeter Brune   PetscReal           fnorm;
930742fe5e2SPeter Brune   SNESConvergedReason reason;
93139bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
93222c1e704SPeter Brune   Vec                 Y, FPC;
933fde0ff24SPeter Brune   PetscBool           lssuccess;
934fde0ff24SPeter Brune   PetscInt            k;
935421d9b32SPeter Brune   PetscFunctionBegin;
936fde0ff24SPeter Brune   Y = snes->work[3];
937d1adcc6fSPeter Brune   if (fas->downsmooth) {
938d1adcc6fSPeter Brune     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
939742fe5e2SPeter Brune     /* check convergence reason for the smoother */
940742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
941742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
942742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
943742fe5e2SPeter Brune       PetscFunctionReturn(0);
944742fe5e2SPeter Brune     }
9454b32a720SPeter Brune     ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
9464b32a720SPeter Brune     ierr = VecCopy(FPC, F);CHKERRQ(ierr);
947fde0ff24SPeter Brune   } else {
948fde0ff24SPeter Brune     /* basic richardson smoother */
949fde0ff24SPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
950794bee33SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
951794bee33SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
952fde0ff24SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
9536188f407SPeter Brune       ierr = PetscLineSearchApply(fas->linesearch_smooth,X,F,&fnorm,Y);CHKERRQ(ierr);
9546188f407SPeter Brune       ierr = PetscLineSearchGetSuccess(fas->linesearch_smooth, &lssuccess);CHKERRQ(ierr);
955fde0ff24SPeter Brune       if (!lssuccess) {
956fde0ff24SPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
957fde0ff24SPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
9582f7ea302SPeter Brune           PetscFunctionReturn(0);
9592f7ea302SPeter Brune         }
960fe6f9142SPeter Brune       }
961fde0ff24SPeter Brune     }
962fde0ff24SPeter Brune   }
96339bd7f45SPeter Brune   PetscFunctionReturn(0);
96439bd7f45SPeter Brune }
96539bd7f45SPeter Brune 
96639bd7f45SPeter Brune 
96739bd7f45SPeter Brune #undef __FUNCT__
96839bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth"
96939bd7f45SPeter Brune /*
97007144faaSPeter Brune Defines the action of the upsmoother
97139bd7f45SPeter Brune  */
97239bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
97339bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
97422c1e704SPeter Brune   PetscReal           fnorm;
97539bd7f45SPeter Brune   SNESConvergedReason reason;
97639bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
97722c1e704SPeter Brune   Vec                 Y, FPC;
978fde0ff24SPeter Brune   PetscBool           lssuccess;
979fde0ff24SPeter Brune   PetscInt            k;
98039bd7f45SPeter Brune   PetscFunctionBegin;
981fde0ff24SPeter Brune   Y = snes->work[3];
98239bd7f45SPeter Brune   if (fas->upsmooth) {
983fde0ff24SPeter Brune     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
98439bd7f45SPeter Brune     /* check convergence reason for the smoother */
985fde0ff24SPeter Brune     ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr);
98639bd7f45SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
98739bd7f45SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
98839bd7f45SPeter Brune       PetscFunctionReturn(0);
98939bd7f45SPeter Brune     }
9904b32a720SPeter Brune     ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
9914b32a720SPeter Brune     ierr = VecCopy(FPC, F);CHKERRQ(ierr);
992fde0ff24SPeter Brune   } else {
993fde0ff24SPeter Brune     /* basic richardson smoother */
994fde0ff24SPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
99539bd7f45SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
99639bd7f45SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
997fde0ff24SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
9986188f407SPeter Brune       ierr = PetscLineSearchApply(fas->linesearch_smooth,X,F,&fnorm,Y);CHKERRQ(ierr);
9996188f407SPeter Brune       ierr = PetscLineSearchGetSuccess(fas->linesearch_smooth, &lssuccess);CHKERRQ(ierr);
1000fde0ff24SPeter Brune       if (!lssuccess) {
1001fde0ff24SPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
1002fde0ff24SPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
100339bd7f45SPeter Brune           PetscFunctionReturn(0);
100439bd7f45SPeter Brune         }
100539bd7f45SPeter Brune       }
1006fde0ff24SPeter Brune     }
1007fde0ff24SPeter Brune   }
100839bd7f45SPeter Brune   PetscFunctionReturn(0);
100939bd7f45SPeter Brune }
101039bd7f45SPeter Brune 
101139bd7f45SPeter Brune #undef __FUNCT__
1012938e4a01SJed Brown #define __FUNCT__ "SNESFASCreateCoarseVec"
1013938e4a01SJed Brown /*@
1014938e4a01SJed Brown    SNESFASCreateCoarseVec - create Vec corresponding to a state vector on one level coarser than current level
1015938e4a01SJed Brown 
1016938e4a01SJed Brown    Collective
1017938e4a01SJed Brown 
1018938e4a01SJed Brown    Input Arguments:
1019938e4a01SJed Brown .  snes - SNESFAS
1020938e4a01SJed Brown 
1021938e4a01SJed Brown    Output Arguments:
1022938e4a01SJed Brown .  Xcoarse - vector on level one coarser than snes
1023938e4a01SJed Brown 
1024938e4a01SJed Brown    Level: developer
1025938e4a01SJed Brown 
1026938e4a01SJed Brown .seealso: SNESFASSetRestriction(), SNESFASRestrict()
1027938e4a01SJed Brown @*/
1028938e4a01SJed Brown PetscErrorCode SNESFASCreateCoarseVec(SNES snes,Vec *Xcoarse)
1029938e4a01SJed Brown {
1030938e4a01SJed Brown   PetscErrorCode ierr;
1031938e4a01SJed Brown   SNES_FAS       *fas = (SNES_FAS*)snes->data;
1032938e4a01SJed Brown 
1033938e4a01SJed Brown   PetscFunctionBegin;
1034938e4a01SJed Brown   if (fas->rscale) {ierr = VecDuplicate(fas->rscale,Xcoarse);CHKERRQ(ierr);}
1035938e4a01SJed Brown   else if (fas->inject) {ierr = MatGetVecs(fas->inject,Xcoarse,PETSC_NULL);CHKERRQ(ierr);}
1036938e4a01SJed Brown   else SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_ARG_WRONGSTATE,"Must set restriction or injection");CHKERRQ(ierr);
1037938e4a01SJed Brown   PetscFunctionReturn(0);
1038938e4a01SJed Brown }
1039938e4a01SJed Brown 
1040e9923e8dSJed Brown #undef __FUNCT__
1041e9923e8dSJed Brown #define __FUNCT__ "SNESFASRestrict"
1042e9923e8dSJed Brown /*@
1043e9923e8dSJed Brown    SNESFASRestrict - restrict a Vec to the next coarser level
1044e9923e8dSJed Brown 
1045e9923e8dSJed Brown    Collective
1046e9923e8dSJed Brown 
1047e9923e8dSJed Brown    Input Arguments:
1048e9923e8dSJed Brown +  fine - SNES from which to restrict
1049e9923e8dSJed Brown -  Xfine - vector to restrict
1050e9923e8dSJed Brown 
1051e9923e8dSJed Brown    Output Arguments:
1052e9923e8dSJed Brown .  Xcoarse - result of restriction
1053e9923e8dSJed Brown 
1054e9923e8dSJed Brown    Level: developer
1055e9923e8dSJed Brown 
1056e9923e8dSJed Brown .seealso: SNESFASSetRestriction(), SNESFASSetInjection()
1057e9923e8dSJed Brown @*/
1058e9923e8dSJed Brown PetscErrorCode SNESFASRestrict(SNES fine,Vec Xfine,Vec Xcoarse)
1059e9923e8dSJed Brown {
1060e9923e8dSJed Brown   PetscErrorCode ierr;
1061e9923e8dSJed Brown   SNES_FAS       *fas = (SNES_FAS*)fine->data;
1062e9923e8dSJed Brown 
1063e9923e8dSJed Brown   PetscFunctionBegin;
1064e9923e8dSJed Brown   PetscValidHeaderSpecific(fine,SNES_CLASSID,1);
1065e9923e8dSJed Brown   PetscValidHeaderSpecific(Xfine,VEC_CLASSID,2);
1066e9923e8dSJed Brown   PetscValidHeaderSpecific(Xcoarse,VEC_CLASSID,3);
1067e9923e8dSJed Brown   if (fas->inject) {
1068e9923e8dSJed Brown     ierr = MatRestrict(fas->inject,Xfine,Xcoarse);CHKERRQ(ierr);
1069e9923e8dSJed Brown   } else {
1070e9923e8dSJed Brown     ierr = MatRestrict(fas->restrct,Xfine,Xcoarse);CHKERRQ(ierr);
1071e9923e8dSJed Brown     ierr = VecPointwiseMult(Xcoarse,fas->rscale,Xcoarse);CHKERRQ(ierr);
1072e9923e8dSJed Brown   }
1073e9923e8dSJed Brown   PetscFunctionReturn(0);
1074e9923e8dSJed Brown }
1075e9923e8dSJed Brown 
1076e9923e8dSJed Brown #undef __FUNCT__
107739bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection"
107839bd7f45SPeter Brune /*
107939bd7f45SPeter Brune 
108039bd7f45SPeter Brune Performs the FAS coarse correction as:
108139bd7f45SPeter Brune 
108239bd7f45SPeter Brune fine problem: F(x) = 0
108339bd7f45SPeter Brune coarse problem: F^c(x) = b^c
108439bd7f45SPeter Brune 
108539bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
108639bd7f45SPeter Brune 
108739bd7f45SPeter Brune with correction:
108839bd7f45SPeter Brune 
108939bd7f45SPeter Brune 
109039bd7f45SPeter Brune 
109139bd7f45SPeter Brune  */
109239a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
109339bd7f45SPeter Brune   PetscErrorCode      ierr;
109439bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
109539bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
109639bd7f45SPeter Brune   SNESConvergedReason reason;
109739bd7f45SPeter Brune   PetscFunctionBegin;
1098fa9694d7SPeter Brune   if (fas->next) {
1099c90fad12SPeter Brune     X_c  = fas->next->vec_sol;
1100293a7e31SPeter Brune     Xo_c = fas->next->work[0];
1101c90fad12SPeter Brune     F_c  = fas->next->vec_func;
1102742fe5e2SPeter Brune     B_c  = fas->next->vec_rhs;
1103efe1f98aSPeter Brune 
1104938e4a01SJed Brown     ierr = SNESFASRestrict(snes,X,Xo_c);CHKERRQ(ierr);
1105293a7e31SPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1106293a7e31SPeter Brune 
1107293a7e31SPeter Brune     /* restrict the defect */
1108293a7e31SPeter Brune     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1109293a7e31SPeter Brune 
1110c90fad12SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1111ee78dd50SPeter Brune     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1112c90fad12SPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1113742fe5e2SPeter Brune 
1114293a7e31SPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1115c90fad12SPeter Brune 
1116ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
1117ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1118c90fad12SPeter Brune 
1119c90fad12SPeter Brune     /* recurse to the next level */
1120f5a6d4f9SBarry Smith     fas->next->vec_rhs = B_c;
1121742fe5e2SPeter Brune     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
1122742fe5e2SPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1123742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1124742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1125742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1126742fe5e2SPeter Brune       PetscFunctionReturn(0);
1127742fe5e2SPeter Brune     }
1128742fe5e2SPeter Brune     /* fas->next->vec_rhs = PETSC_NULL; */
1129ee78dd50SPeter Brune 
1130fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
1131fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
113239bd7f45SPeter Brune     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr);
1133293a7e31SPeter Brune   }
113439bd7f45SPeter Brune   PetscFunctionReturn(0);
113539bd7f45SPeter Brune }
113639bd7f45SPeter Brune 
113739bd7f45SPeter Brune #undef __FUNCT__
113839bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive"
113939bd7f45SPeter Brune /*
114039bd7f45SPeter Brune 
114139bd7f45SPeter Brune The additive cycle looks like:
114239bd7f45SPeter Brune 
114307144faaSPeter Brune xhat = x
114407144faaSPeter Brune xhat = dS(x, b)
114507144faaSPeter Brune x = coarsecorrection(xhat, b_d)
114607144faaSPeter Brune x = x + nu*(xhat - x);
114739bd7f45SPeter Brune (optional) x = uS(x, b)
114839bd7f45SPeter Brune 
114939bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
115039bd7f45SPeter Brune 
115139bd7f45SPeter Brune  */
115239bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
115307144faaSPeter Brune   Vec                 F, B, Xhat;
115422c1e704SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
115539bd7f45SPeter Brune   PetscErrorCode      ierr;
115607144faaSPeter Brune   SNES_FAS *          fas = (SNES_FAS *)snes->data;
115707144faaSPeter Brune   SNESConvergedReason reason;
115822c1e704SPeter Brune   PetscReal           xnorm, fnorm, ynorm;
115922c1e704SPeter Brune   PetscBool           lssuccess;
116039bd7f45SPeter Brune   PetscFunctionBegin;
116139bd7f45SPeter Brune 
116239bd7f45SPeter Brune   F = snes->vec_func;
116339bd7f45SPeter Brune   B = snes->vec_rhs;
1164fde0ff24SPeter Brune   Xhat = snes->work[3];
116507144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
116607144faaSPeter Brune   /* recurse first */
116707144faaSPeter Brune   if (fas->next) {
116807144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
116939bd7f45SPeter Brune 
117007144faaSPeter Brune     X_c  = fas->next->vec_sol;
117107144faaSPeter Brune     Xo_c = fas->next->work[0];
117207144faaSPeter Brune     F_c  = fas->next->vec_func;
117307144faaSPeter Brune     B_c  = fas->next->vec_rhs;
117439bd7f45SPeter Brune 
1175938e4a01SJed Brown     ierr = SNESFASRestrict(snes,Xhat,Xo_c);CHKERRQ(ierr);
117607144faaSPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
117707144faaSPeter Brune 
117807144faaSPeter Brune     /* restrict the defect */
117907144faaSPeter Brune     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
118007144faaSPeter Brune 
118107144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
118207144faaSPeter Brune     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
118307144faaSPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
118407144faaSPeter Brune 
118507144faaSPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
118607144faaSPeter Brune 
118707144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
118807144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
118907144faaSPeter Brune 
119007144faaSPeter Brune     /* recurse */
119107144faaSPeter Brune     fas->next->vec_rhs = B_c;
119207144faaSPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
119307144faaSPeter Brune 
119407144faaSPeter Brune     /* smooth on this level */
119507144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
119607144faaSPeter Brune 
119707144faaSPeter Brune     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
119807144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
119907144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
120007144faaSPeter Brune       PetscFunctionReturn(0);
120107144faaSPeter Brune     }
120207144faaSPeter Brune 
120307144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
1204c68acad4SPeter Brune     ierr = VecAYPX(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1205ddebd997SPeter Brune     ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr);
120607144faaSPeter Brune 
1207ddebd997SPeter Brune     /* additive correction of the coarse direction*/
1208ddebd997SPeter Brune     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1209ddebd997SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
12109e764e56SPeter Brune     ierr = PetscLineSearchApply(snes->linesearch, X, F, &fnorm, Xhat);CHKERRQ(ierr);
12119e764e56SPeter Brune     ierr = PetscLineSearchGetSuccess(snes->linesearch, &lssuccess);CHKERRQ(ierr);
12129e764e56SPeter Brune     if (!lssuccess) {
12139e764e56SPeter Brune       if (++snes->numFailures >= snes->maxFailures) {
12149e764e56SPeter Brune         snes->reason = SNES_DIVERGED_LINE_SEARCH;
12159e764e56SPeter Brune         PetscFunctionReturn(0);
12169e764e56SPeter Brune       }
12179e764e56SPeter Brune     }
12189e764e56SPeter Brune     ierr = PetscLineSearchGetNorms(snes->linesearch, &xnorm, &fnorm, &ynorm);CHKERRQ(ierr);
121907144faaSPeter Brune   } else {
122007144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
122107144faaSPeter Brune   }
122239bd7f45SPeter Brune   PetscFunctionReturn(0);
122339bd7f45SPeter Brune }
122439bd7f45SPeter Brune 
122539bd7f45SPeter Brune #undef __FUNCT__
122639bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative"
122739bd7f45SPeter Brune /*
122839bd7f45SPeter Brune 
122939bd7f45SPeter Brune Defines the FAS cycle as:
123039bd7f45SPeter Brune 
123139bd7f45SPeter Brune fine problem: F(x) = 0
123239bd7f45SPeter Brune coarse problem: F^c(x) = b^c
123339bd7f45SPeter Brune 
123439bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
123539bd7f45SPeter Brune 
123639bd7f45SPeter Brune correction:
123739bd7f45SPeter Brune 
123839bd7f45SPeter Brune x = x + I(x^c - Rx)
123939bd7f45SPeter Brune 
124039bd7f45SPeter Brune  */
124139bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
124239bd7f45SPeter Brune 
124339bd7f45SPeter Brune   PetscErrorCode      ierr;
124439bd7f45SPeter Brune   Vec                 F,B;
124539bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
124639bd7f45SPeter Brune 
124739bd7f45SPeter Brune   PetscFunctionBegin;
124839bd7f45SPeter Brune   F = snes->vec_func;
124939bd7f45SPeter Brune   B = snes->vec_rhs;
125039bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
125139bd7f45SPeter Brune   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
125239bd7f45SPeter Brune 
125339a4b097SPeter Brune   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
125439bd7f45SPeter Brune 
1255c90fad12SPeter Brune   if (fas->level != 0) {
125639bd7f45SPeter Brune     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
1257fe6f9142SPeter Brune   }
1258fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1259fa9694d7SPeter Brune 
1260fa9694d7SPeter Brune   PetscFunctionReturn(0);
1261421d9b32SPeter Brune }
1262421d9b32SPeter Brune 
1263421d9b32SPeter Brune #undef __FUNCT__
1264421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
1265421d9b32SPeter Brune 
1266421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
1267421d9b32SPeter Brune {
1268fa9694d7SPeter Brune   PetscErrorCode ierr;
1269fe6f9142SPeter Brune   PetscInt       i, maxits;
1270ddb5aff1SPeter Brune   Vec            X, F;
1271fe6f9142SPeter Brune   PetscReal      fnorm;
1272b17ce1afSJed Brown   SNES_FAS       *fas = (SNES_FAS *)snes->data,*ffas;
1273b17ce1afSJed Brown   DM             dm;
1274b17ce1afSJed Brown 
1275421d9b32SPeter Brune   PetscFunctionBegin;
1276fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
1277fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1278fa9694d7SPeter Brune   X = snes->vec_sol;
1279f5a6d4f9SBarry Smith   F = snes->vec_func;
1280293a7e31SPeter Brune 
1281293a7e31SPeter Brune   /*norm setup */
1282fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1283fe6f9142SPeter Brune   snes->iter = 0;
1284fe6f9142SPeter Brune   snes->norm = 0.;
1285fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1286fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1287fe6f9142SPeter Brune   if (snes->domainerror) {
1288fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1289fe6f9142SPeter Brune     PetscFunctionReturn(0);
1290fe6f9142SPeter Brune   }
1291fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1292fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1293fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1294fe6f9142SPeter Brune   snes->norm = fnorm;
1295fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1296fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1297fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1298fe6f9142SPeter Brune 
1299fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
1300fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
1301fe6f9142SPeter Brune   /* test convergence */
1302fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1303fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
1304b17ce1afSJed Brown 
1305b17ce1afSJed Brown   ierr = SNESGetDM(snes,&dm);CHKERRQ(ierr);
1306b17ce1afSJed Brown   for (ffas=fas; ffas->next; ffas=(SNES_FAS*)ffas->next->data) {
1307b17ce1afSJed Brown     DM dmcoarse;
1308b17ce1afSJed Brown     ierr = SNESGetDM(ffas->next,&dmcoarse);CHKERRQ(ierr);
1309b17ce1afSJed Brown     ierr = DMRestrict(dm,ffas->restrct,ffas->rscale,ffas->inject,dmcoarse);CHKERRQ(ierr);
1310b17ce1afSJed Brown     dm = dmcoarse;
1311b17ce1afSJed Brown   }
1312b17ce1afSJed Brown 
1313fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
1314fe6f9142SPeter Brune     /* Call general purpose update function */
1315646217ecSPeter Brune 
1316fe6f9142SPeter Brune     if (snes->ops->update) {
1317fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1318fe6f9142SPeter Brune     }
131907144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
132039bd7f45SPeter Brune       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
132107144faaSPeter Brune     } else {
132207144faaSPeter Brune       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
132307144faaSPeter Brune     }
1324742fe5e2SPeter Brune 
1325742fe5e2SPeter Brune     /* check for FAS cycle divergence */
1326742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
1327742fe5e2SPeter Brune       PetscFunctionReturn(0);
1328742fe5e2SPeter Brune     }
1329c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1330c90fad12SPeter Brune     /* Monitor convergence */
1331c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1332c90fad12SPeter Brune     snes->iter = i+1;
1333c90fad12SPeter Brune     snes->norm = fnorm;
1334c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1335c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
1336c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1337c90fad12SPeter Brune     /* Test for convergence */
1338c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1339c90fad12SPeter Brune     if (snes->reason) break;
1340fe6f9142SPeter Brune   }
1341fe6f9142SPeter Brune   if (i == maxits) {
1342fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1343fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1344fe6f9142SPeter Brune   }
1345421d9b32SPeter Brune   PetscFunctionReturn(0);
1346421d9b32SPeter Brune }
1347