xref: /petsc/src/snes/impls/fas/fas.c (revision 07144faae8d9ae592719efdbb82dce8f06715f56)
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 /*MC
7421d9b32SPeter Brune Full Approximation Scheme nonlinear multigrid solver.
8421d9b32SPeter Brune 
9421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
10421d9b32SPeter Brune 
11421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
12421d9b32SPeter Brune M*/
13421d9b32SPeter Brune 
14421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
15421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
16421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
17421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
18421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
19421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
206273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *);
21421d9b32SPeter Brune 
22421d9b32SPeter Brune EXTERN_C_BEGIN
23421d9b32SPeter Brune 
24421d9b32SPeter Brune #undef __FUNCT__
25421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
26421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes)
27421d9b32SPeter Brune {
28421d9b32SPeter Brune   SNES_FAS * fas;
29421d9b32SPeter Brune   PetscErrorCode ierr;
30421d9b32SPeter Brune 
31421d9b32SPeter Brune   PetscFunctionBegin;
32421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
33421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
34421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
35421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
36421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
37421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
38421d9b32SPeter Brune 
39ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
40ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
41ed020824SBarry Smith 
42421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
43421d9b32SPeter Brune   snes->data                  = (void*) fas;
44421d9b32SPeter Brune   fas->level                  = 0;
45293a7e31SPeter Brune   fas->levels                 = 1;
46ee78dd50SPeter Brune   fas->n_cycles               = 1;
47ee78dd50SPeter Brune   fas->max_up_it              = 1;
48ee78dd50SPeter Brune   fas->max_down_it            = 1;
49ee78dd50SPeter Brune   fas->upsmooth               = PETSC_NULL;
50ee78dd50SPeter Brune   fas->downsmooth             = PETSC_NULL;
51421d9b32SPeter Brune   fas->next                   = PETSC_NULL;
526273346dSPeter Brune   fas->previous               = PETSC_NULL;
53421d9b32SPeter Brune   fas->interpolate            = PETSC_NULL;
54421d9b32SPeter Brune   fas->restrct                = PETSC_NULL;
55efe1f98aSPeter Brune   fas->inject                 = PETSC_NULL;
56cc05f883SPeter Brune   fas->monitor                = PETSC_NULL;
57cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
58421d9b32SPeter Brune   PetscFunctionReturn(0);
59421d9b32SPeter Brune }
60421d9b32SPeter Brune EXTERN_C_END
61421d9b32SPeter Brune 
62421d9b32SPeter Brune #undef __FUNCT__
63421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
64c79ef259SPeter Brune /*@
65722262beSPeter Brune    SNESFASGetLevels - Gets the number of levels in a FAS.
66c79ef259SPeter Brune 
67c79ef259SPeter Brune    Input Parameter:
68c79ef259SPeter Brune .  snes - the preconditioner context
69c79ef259SPeter Brune 
70c79ef259SPeter Brune    Output parameter:
71c79ef259SPeter Brune .  levels - the number of levels
72c79ef259SPeter Brune 
73c79ef259SPeter Brune    Level: advanced
74c79ef259SPeter Brune 
75c79ef259SPeter Brune .keywords: MG, get, levels, multigrid
76c79ef259SPeter Brune 
77c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels()
78c79ef259SPeter Brune @*/
79421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
80421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
81421d9b32SPeter Brune   PetscFunctionBegin;
82ee1fd11aSPeter Brune   *levels = fas->levels;
83421d9b32SPeter Brune   PetscFunctionReturn(0);
84421d9b32SPeter Brune }
85421d9b32SPeter Brune 
86421d9b32SPeter Brune #undef __FUNCT__
87646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles"
88c79ef259SPeter Brune /*@
89c79ef259SPeter Brune    SNESFASSetCycles - Sets the type cycles to use.  Use SNESFASSetCyclesOnLevel() for more
90c79ef259SPeter Brune    complicated cycling.
91c79ef259SPeter Brune 
92c79ef259SPeter Brune    Logically Collective on SNES
93c79ef259SPeter Brune 
94c79ef259SPeter Brune    Input Parameters:
95c79ef259SPeter Brune +  snes   - the multigrid context
96c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
97c79ef259SPeter Brune 
98c79ef259SPeter Brune    Options Database Key:
99c79ef259SPeter Brune $  -snes_fas_cycles 1 or 2
100c79ef259SPeter Brune 
101c79ef259SPeter Brune    Level: advanced
102c79ef259SPeter Brune 
103c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
104c79ef259SPeter Brune 
105c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel()
106c79ef259SPeter Brune @*/
107646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
108646217ecSPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
109646217ecSPeter Brune   PetscErrorCode ierr;
110646217ecSPeter Brune   PetscFunctionBegin;
111646217ecSPeter Brune   fas->n_cycles = cycles;
112646217ecSPeter Brune   if (fas->next) {
113646217ecSPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
114646217ecSPeter Brune   }
115646217ecSPeter Brune   PetscFunctionReturn(0);
116646217ecSPeter Brune }
117646217ecSPeter Brune 
118eff52c0eSPeter Brune #undef __FUNCT__
119c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel"
120c79ef259SPeter Brune /*@
121722262beSPeter Brune    SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level.
122c79ef259SPeter Brune 
123c79ef259SPeter Brune    Logically Collective on SNES
124c79ef259SPeter Brune 
125c79ef259SPeter Brune    Input Parameters:
126c79ef259SPeter Brune +  snes   - the multigrid context
127c79ef259SPeter Brune .  level  - the level to set the number of cycles on
128c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
129c79ef259SPeter Brune 
130c79ef259SPeter Brune    Level: advanced
131c79ef259SPeter Brune 
132c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
133c79ef259SPeter Brune 
134c79ef259SPeter Brune .seealso: SNESFASSetCycles()
135c79ef259SPeter Brune @*/
136c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) {
137c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
138c79ef259SPeter Brune   PetscInt top_level = fas->level,i;
139c79ef259SPeter Brune 
140c79ef259SPeter Brune   PetscFunctionBegin;
141c79ef259SPeter Brune   if (level > top_level)
142c79ef259SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
143c79ef259SPeter Brune   /* get to the correct level */
144c79ef259SPeter Brune   for (i = fas->level; i > level; i--) {
145c79ef259SPeter Brune     fas = (SNES_FAS *)fas->next->data;
146c79ef259SPeter Brune   }
147c79ef259SPeter Brune   if (fas->level != level)
148c79ef259SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
149c79ef259SPeter Brune   fas->n_cycles = cycles;
150c79ef259SPeter Brune   PetscFunctionReturn(0);
151c79ef259SPeter Brune }
152c79ef259SPeter Brune 
153c79ef259SPeter Brune #undef __FUNCT__
154eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS"
155aeed3662SMatthew G Knepley /*@C
156c79ef259SPeter Brune    SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used.
157c79ef259SPeter Brune    Use SNESFASSetGSOnLevel() for more complicated staging of smoothers
158c79ef259SPeter Brune    and nonlinear preconditioners.
159c79ef259SPeter Brune 
160c79ef259SPeter Brune    Logically Collective on SNES
161c79ef259SPeter Brune 
162c79ef259SPeter Brune    Input Parameters:
163c79ef259SPeter Brune +  snes    - the multigrid context
164c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
165c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
166c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
167c79ef259SPeter Brune 
168c79ef259SPeter Brune    Level: advanced
169c79ef259SPeter Brune 
170c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid
171c79ef259SPeter Brune 
172c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel()
173c79ef259SPeter Brune @*/
174eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
175eff52c0eSPeter Brune   PetscErrorCode ierr = 0;
176eff52c0eSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
177eff52c0eSPeter Brune   PetscFunctionBegin;
178eff52c0eSPeter Brune 
179eff52c0eSPeter Brune   /* use or don't use it according to user wishes*/
180d28543b3SPeter Brune   snes->usegs = use_gs;
181eff52c0eSPeter Brune   if (gsfunc) {
182eff52c0eSPeter Brune     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
183eff52c0eSPeter Brune     /* push the provided GS up the tree */
184eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
185eff52c0eSPeter Brune   } else if (snes->ops->computegs) {
186eff52c0eSPeter Brune     /* assume that the user has set the GS solver at this level */
187eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
188eff52c0eSPeter Brune   } else if (use_gs) {
189eff52c0eSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level);
190eff52c0eSPeter Brune   }
191eff52c0eSPeter Brune   PetscFunctionReturn(0);
192eff52c0eSPeter Brune }
193eff52c0eSPeter Brune 
194eff52c0eSPeter Brune #undef __FUNCT__
195eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel"
196aeed3662SMatthew G Knepley /*@C
197c79ef259SPeter Brune    SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level.
198c79ef259SPeter Brune 
199c79ef259SPeter Brune    Logically Collective on SNES
200c79ef259SPeter Brune 
201c79ef259SPeter Brune    Input Parameters:
202c79ef259SPeter Brune +  snes    - the multigrid context
203c79ef259SPeter Brune .  level   - the level to set the nonlinear smoother on
204c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
205c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
206c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
207c79ef259SPeter Brune 
208c79ef259SPeter Brune    Level: advanced
209c79ef259SPeter Brune 
210c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
211c79ef259SPeter Brune 
212c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS()
213c79ef259SPeter Brune @*/
214eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
215eff52c0eSPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
216eff52c0eSPeter Brune   PetscErrorCode ierr;
217eff52c0eSPeter Brune   PetscInt       top_level = fas->level,i;
218eff52c0eSPeter Brune   SNES           cur_snes = snes;
219eff52c0eSPeter Brune   PetscFunctionBegin;
220eff52c0eSPeter Brune   if (level > top_level)
221eff52c0eSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
222eff52c0eSPeter Brune   /* get to the correct level */
223eff52c0eSPeter Brune   for (i = fas->level; i > level; i--) {
224eff52c0eSPeter Brune     fas = (SNES_FAS *)fas->next->data;
225eff52c0eSPeter Brune     cur_snes = fas->next;
226eff52c0eSPeter Brune   }
227eff52c0eSPeter Brune   if (fas->level != level)
228eff52c0eSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
229d28543b3SPeter Brune   snes->usegs = use_gs;
230eff52c0eSPeter Brune   if (gsfunc) {
2316273346dSPeter Brune     ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr);
232eff52c0eSPeter Brune   }
233eff52c0eSPeter Brune   PetscFunctionReturn(0);
234eff52c0eSPeter Brune }
235eff52c0eSPeter Brune 
236646217ecSPeter Brune #undef __FUNCT__
237421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
238c79ef259SPeter Brune /*@
239c79ef259SPeter Brune    SNESFASGetSNES - Gets the SNES corresponding to a particular
240c79ef259SPeter Brune    level of the FAS hierarchy.
241c79ef259SPeter Brune 
242c79ef259SPeter Brune    Input Parameters:
243c79ef259SPeter Brune +  snes    - the multigrid context
244c79ef259SPeter Brune    level   - the level to get
245c79ef259SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
246c79ef259SPeter Brune 
247c79ef259SPeter Brune    Level: advanced
248c79ef259SPeter Brune 
249c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
250c79ef259SPeter Brune 
251c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
252c79ef259SPeter Brune @*/
253421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
254421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
255421d9b32SPeter Brune   PetscInt levels = fas->level;
256421d9b32SPeter Brune   PetscInt i;
257421d9b32SPeter Brune   PetscFunctionBegin;
258421d9b32SPeter Brune   *lsnes = snes;
259421d9b32SPeter Brune   if (fas->level < level) {
260421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
261421d9b32SPeter Brune   }
262421d9b32SPeter Brune   if (level > levels - 1) {
263421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
264421d9b32SPeter Brune   }
265421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
266421d9b32SPeter Brune     *lsnes = fas->next;
267421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
268421d9b32SPeter Brune   }
269421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
270421d9b32SPeter Brune   PetscFunctionReturn(0);
271421d9b32SPeter Brune }
272421d9b32SPeter Brune 
273421d9b32SPeter Brune #undef __FUNCT__
27407144faaSPeter Brune #define __FUNCT__ "SNESFASSetType"
27507144faaSPeter Brune /*@
27607144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS.
27707144faaSPeter Brune e
27807144faaSPeter Brune 
27907144faaSPeter Brune 
28007144faaSPeter Brune @*/
28107144faaSPeter Brune PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
28207144faaSPeter Brune {
28307144faaSPeter Brune   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
28407144faaSPeter Brune 
28507144faaSPeter Brune   PetscFunctionBegin;
28607144faaSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
28707144faaSPeter Brune   PetscValidLogicalCollectiveEnum(snes,fastype,2);
28807144faaSPeter Brune   fas->fastype = fastype;
28907144faaSPeter Brune   PetscFunctionReturn(0);
29007144faaSPeter Brune }
29107144faaSPeter Brune 
29207144faaSPeter Brune 
29307144faaSPeter Brune 
29407144faaSPeter Brune #undef __FUNCT__
295421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
296aeed3662SMatthew G Knepley /*@C
297c79ef259SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
298c79ef259SPeter Brune    Must be called before any other FAS routine.
299c79ef259SPeter Brune 
300c79ef259SPeter Brune    Input Parameters:
301c79ef259SPeter Brune +  snes   - the snes context
302c79ef259SPeter Brune .  levels - the number of levels
303c79ef259SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
304c79ef259SPeter Brune             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
305c79ef259SPeter Brune             Fortran.
306c79ef259SPeter Brune 
307c79ef259SPeter Brune    Level: intermediate
308c79ef259SPeter Brune 
309c79ef259SPeter Brune    Notes:
310c79ef259SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
311c79ef259SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
312c79ef259SPeter Brune 
313c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid
314c79ef259SPeter Brune 
315c79ef259SPeter Brune .seealso: SNESFASGetLevels()
316c79ef259SPeter Brune @*/
317421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
318421d9b32SPeter Brune   PetscErrorCode ierr;
319421d9b32SPeter Brune   PetscInt i;
320421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
3216273346dSPeter Brune   SNES prevsnes;
322421d9b32SPeter Brune   MPI_Comm comm;
323421d9b32SPeter Brune   PetscFunctionBegin;
324ee1fd11aSPeter Brune   comm = ((PetscObject)snes)->comm;
325ee1fd11aSPeter Brune   if (levels == fas->levels) {
326ee1fd11aSPeter Brune     if (!comms)
327ee1fd11aSPeter Brune       PetscFunctionReturn(0);
328ee1fd11aSPeter Brune   }
329421d9b32SPeter Brune   /* user has changed the number of levels; reset */
330421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
331421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
332421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
333ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
3346273346dSPeter Brune   fas->previous = PETSC_NULL;
3356273346dSPeter Brune   prevsnes = snes;
336421d9b32SPeter Brune   /* setup the finest level */
337421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
338421d9b32SPeter Brune     if (comms) comm = comms[i];
339421d9b32SPeter Brune     fas->level = i;
340421d9b32SPeter Brune     fas->levels = levels;
341ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
342421d9b32SPeter Brune     if (i > 0) {
343421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
3444a6b3fb8SBarry Smith       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
345421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
346794bee33SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
3476273346dSPeter Brune       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
3486273346dSPeter Brune       prevsnes = fas->next;
3496273346dSPeter Brune       fas = (SNES_FAS *)prevsnes->data;
350421d9b32SPeter Brune     }
351421d9b32SPeter Brune   }
352421d9b32SPeter Brune   PetscFunctionReturn(0);
353421d9b32SPeter Brune }
354421d9b32SPeter Brune 
355421d9b32SPeter Brune #undef __FUNCT__
356c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp"
357c79ef259SPeter Brune /*@
358c79ef259SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
359c79ef259SPeter Brune    use on all levels.
360c79ef259SPeter Brune 
361c79ef259SPeter Brune    Logically Collective on PC
362c79ef259SPeter Brune 
363c79ef259SPeter Brune    Input Parameters:
364c79ef259SPeter Brune +  snes - the multigrid context
365c79ef259SPeter Brune -  n    - the number of smoothing steps
366c79ef259SPeter Brune 
367c79ef259SPeter Brune    Options Database Key:
368d28543b3SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
369c79ef259SPeter Brune 
370c79ef259SPeter Brune    Level: advanced
371c79ef259SPeter Brune 
372c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
373c79ef259SPeter Brune 
374c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
375c79ef259SPeter Brune @*/
376c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
377c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
378c79ef259SPeter Brune   PetscErrorCode ierr = 0;
379c79ef259SPeter Brune   PetscFunctionBegin;
380c79ef259SPeter Brune   fas->max_up_it = n;
381c79ef259SPeter Brune   if (fas->next) {
382c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
383c79ef259SPeter Brune   }
384c79ef259SPeter Brune   PetscFunctionReturn(0);
385c79ef259SPeter Brune }
386c79ef259SPeter Brune 
387c79ef259SPeter Brune #undef __FUNCT__
388c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown"
389c79ef259SPeter Brune /*@
390c79ef259SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
391c79ef259SPeter Brune    use on all levels.
392c79ef259SPeter Brune 
393c79ef259SPeter Brune    Logically Collective on PC
394c79ef259SPeter Brune 
395c79ef259SPeter Brune    Input Parameters:
396c79ef259SPeter Brune +  snes - the multigrid context
397c79ef259SPeter Brune -  n    - the number of smoothing steps
398c79ef259SPeter Brune 
399c79ef259SPeter Brune    Options Database Key:
400d28543b3SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
401c79ef259SPeter Brune 
402c79ef259SPeter Brune    Level: advanced
403c79ef259SPeter Brune 
404c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
405c79ef259SPeter Brune 
406c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
407c79ef259SPeter Brune @*/
408c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
409c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
410c79ef259SPeter Brune   PetscErrorCode ierr = 0;
411c79ef259SPeter Brune   PetscFunctionBegin;
412c79ef259SPeter Brune   fas->max_down_it = n;
413c79ef259SPeter Brune   if (fas->next) {
414c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
415c79ef259SPeter Brune   }
416c79ef259SPeter Brune   PetscFunctionReturn(0);
417c79ef259SPeter Brune }
418c79ef259SPeter Brune 
419c79ef259SPeter Brune #undef __FUNCT__
420421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
421c79ef259SPeter Brune /*@
422c79ef259SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
423c79ef259SPeter Brune    interpolation from l-1 to the lth level
424c79ef259SPeter Brune 
425c79ef259SPeter Brune    Input Parameters:
426c79ef259SPeter Brune +  snes      - the multigrid context
427c79ef259SPeter Brune .  mat       - the interpolation operator
428c79ef259SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
429c79ef259SPeter Brune 
430c79ef259SPeter Brune    Level: advanced
431c79ef259SPeter Brune 
432c79ef259SPeter Brune    Notes:
433c79ef259SPeter Brune           Usually this is the same matrix used also to set the restriction
434c79ef259SPeter Brune     for the same level.
435c79ef259SPeter Brune 
436c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
437c79ef259SPeter Brune     out from the matrix size which one it is.
438c79ef259SPeter Brune 
439c79ef259SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
440c79ef259SPeter Brune 
441c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale()
442c79ef259SPeter Brune @*/
443421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
444421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
445d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
446421d9b32SPeter Brune 
447421d9b32SPeter Brune   PetscFunctionBegin;
448421d9b32SPeter Brune   if (level > top_level)
449421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
450421d9b32SPeter Brune   /* get to the correct level */
451d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
452421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
453421d9b32SPeter Brune   }
454421d9b32SPeter Brune   if (fas->level != level)
455421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
456421d9b32SPeter Brune   fas->interpolate = mat;
457421d9b32SPeter Brune   PetscFunctionReturn(0);
458421d9b32SPeter Brune }
459421d9b32SPeter Brune 
460421d9b32SPeter Brune #undef __FUNCT__
461421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
462c79ef259SPeter Brune /*@
463c79ef259SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
464c79ef259SPeter Brune    from level l to l-1.
465c79ef259SPeter Brune 
466c79ef259SPeter Brune    Input Parameters:
467c79ef259SPeter Brune +  snes  - the multigrid context
468c79ef259SPeter Brune .  mat   - the restriction matrix
469c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
470c79ef259SPeter Brune 
471c79ef259SPeter Brune    Level: advanced
472c79ef259SPeter Brune 
473c79ef259SPeter Brune    Notes:
474c79ef259SPeter Brune           Usually this is the same matrix used also to set the interpolation
475c79ef259SPeter Brune     for the same level.
476c79ef259SPeter Brune 
477c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
478c79ef259SPeter Brune     out from the matrix size which one it is.
479c79ef259SPeter Brune 
480c79ef259SPeter Brune          If you do not set this, the transpose of the Mat set with PCMGSetInterpolation()
481c79ef259SPeter Brune     is used.
482c79ef259SPeter Brune 
483c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
484c79ef259SPeter Brune 
485c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
486c79ef259SPeter Brune @*/
487421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
488421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
489d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
490421d9b32SPeter Brune 
491421d9b32SPeter Brune   PetscFunctionBegin;
492421d9b32SPeter Brune   if (level > top_level)
493421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
494421d9b32SPeter Brune   /* get to the correct level */
495d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
496421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
497421d9b32SPeter Brune   }
498421d9b32SPeter Brune   if (fas->level != level)
499421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
500421d9b32SPeter Brune   fas->restrct = mat;
501421d9b32SPeter Brune   PetscFunctionReturn(0);
502421d9b32SPeter Brune }
503421d9b32SPeter Brune 
504efe1f98aSPeter Brune 
505421d9b32SPeter Brune #undef __FUNCT__
506421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
507c79ef259SPeter Brune /*@
508c79ef259SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
509c79ef259SPeter Brune    operator from level l to l-1.
510c79ef259SPeter Brune 
511c79ef259SPeter Brune    Input Parameters:
512c79ef259SPeter Brune +  snes   - the multigrid context
513c79ef259SPeter Brune .  rscale - the restriction scaling
514c79ef259SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
515c79ef259SPeter Brune 
516c79ef259SPeter Brune    Level: advanced
517c79ef259SPeter Brune 
518c79ef259SPeter Brune    Notes:
519c79ef259SPeter Brune          This is only used in the case that the injection is not set.
520c79ef259SPeter Brune 
521c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
522c79ef259SPeter Brune 
523c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
524c79ef259SPeter Brune @*/
525421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
526421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
527d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
528421d9b32SPeter Brune 
529421d9b32SPeter Brune   PetscFunctionBegin;
530421d9b32SPeter Brune   if (level > top_level)
531421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
532421d9b32SPeter Brune   /* get to the correct level */
533d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
534421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
535421d9b32SPeter Brune   }
536421d9b32SPeter Brune   if (fas->level != level)
537421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
538421d9b32SPeter Brune   fas->rscale = rscale;
539421d9b32SPeter Brune   PetscFunctionReturn(0);
540421d9b32SPeter Brune }
541421d9b32SPeter Brune 
542efe1f98aSPeter Brune 
543efe1f98aSPeter Brune #undef __FUNCT__
544efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection"
545c79ef259SPeter Brune /*@
546c79ef259SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
547c79ef259SPeter Brune    from level l to l-1.
548c79ef259SPeter Brune 
549c79ef259SPeter Brune    Input Parameters:
550c79ef259SPeter Brune +  snes  - the multigrid context
551c79ef259SPeter Brune .  mat   - the restriction matrix
552c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
553c79ef259SPeter Brune 
554c79ef259SPeter Brune    Level: advanced
555c79ef259SPeter Brune 
556c79ef259SPeter Brune    Notes:
557c79ef259SPeter Brune          If you do not set this, the restriction and rscale is used to
558c79ef259SPeter Brune    project the solution instead.
559c79ef259SPeter Brune 
560c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
561c79ef259SPeter Brune 
562c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
563c79ef259SPeter Brune @*/
564efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
565efe1f98aSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
566efe1f98aSPeter Brune   PetscInt top_level = fas->level,i;
567efe1f98aSPeter Brune 
568efe1f98aSPeter Brune   PetscFunctionBegin;
569efe1f98aSPeter Brune   if (level > top_level)
570efe1f98aSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
571efe1f98aSPeter Brune   /* get to the correct level */
572efe1f98aSPeter Brune   for (i = fas->level; i > level; i--) {
573efe1f98aSPeter Brune     fas = (SNES_FAS *)fas->next->data;
574efe1f98aSPeter Brune   }
575efe1f98aSPeter Brune   if (fas->level != level)
576efe1f98aSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
577efe1f98aSPeter Brune   fas->inject = mat;
578efe1f98aSPeter Brune   PetscFunctionReturn(0);
579efe1f98aSPeter Brune }
580efe1f98aSPeter Brune 
581421d9b32SPeter Brune #undef __FUNCT__
582421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
583421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
584421d9b32SPeter Brune {
58577df8cc4SPeter Brune   PetscErrorCode ierr = 0;
586421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
587421d9b32SPeter Brune 
588421d9b32SPeter Brune   PetscFunctionBegin;
589742fe5e2SPeter Brune   if (fas->upsmooth)   ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
590742fe5e2SPeter Brune   if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
591742fe5e2SPeter Brune   if (fas->next)       ierr = SNESReset(fas->next);CHKERRQ(ierr);
592421d9b32SPeter Brune   PetscFunctionReturn(0);
593421d9b32SPeter Brune }
594421d9b32SPeter Brune 
595421d9b32SPeter Brune #undef __FUNCT__
596421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
597421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
598421d9b32SPeter Brune {
599421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
600742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
601421d9b32SPeter Brune 
602421d9b32SPeter Brune   PetscFunctionBegin;
603421d9b32SPeter Brune   /* recursively resets and then destroys */
60451e86f29SPeter Brune   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
60551e86f29SPeter Brune   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
606efe1f98aSPeter Brune   if (fas->inject) {
607efe1f98aSPeter Brune     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
608efe1f98aSPeter Brune   }
60951e86f29SPeter Brune   if (fas->interpolate == fas->restrct) {
61051e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
61151e86f29SPeter Brune     fas->restrct = PETSC_NULL;
61251e86f29SPeter Brune   } else {
61351e86f29SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
61451e86f29SPeter Brune     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
61551e86f29SPeter Brune   }
61651e86f29SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
617421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
618421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
619ee1fd11aSPeter Brune 
620421d9b32SPeter Brune   PetscFunctionReturn(0);
621421d9b32SPeter Brune }
622421d9b32SPeter Brune 
623421d9b32SPeter Brune #undef __FUNCT__
624421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
625421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
626421d9b32SPeter Brune {
62748bfdf8aSPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
628421d9b32SPeter Brune   PetscErrorCode ierr;
629efe1f98aSPeter Brune   VecScatter     injscatter;
630d1adcc6fSPeter Brune   PetscInt       dm_levels;
631d1adcc6fSPeter Brune 
632421d9b32SPeter Brune   PetscFunctionBegin;
633eff52c0eSPeter Brune 
634cc05f883SPeter Brune   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
635d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
636d1adcc6fSPeter Brune     dm_levels++;
637cc05f883SPeter Brune     if (dm_levels > fas->levels) {
638d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
639cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
640d1adcc6fSPeter Brune     }
641d1adcc6fSPeter Brune   }
642d1adcc6fSPeter Brune 
64307144faaSPeter Brune   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
6442f7ea302SPeter Brune     ierr = SNESDefaultGetWork(snes, 1);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
64507144faaSPeter Brune   } else {
64607144faaSPeter Brune     ierr = SNESDefaultGetWork(snes, 2);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
64707144faaSPeter Brune   }
648cc05f883SPeter Brune 
64948bfdf8aSPeter Brune   /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */
65048bfdf8aSPeter Brune   if (fas->upsmooth) {
65148bfdf8aSPeter Brune     ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);
65248bfdf8aSPeter Brune     if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) {
65348bfdf8aSPeter Brune       ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
65448bfdf8aSPeter Brune     }
65548bfdf8aSPeter Brune     if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) {
65648bfdf8aSPeter Brune       ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
65748bfdf8aSPeter Brune     }
65848bfdf8aSPeter Brune    if (snes->ops->computegs && !fas->upsmooth->ops->computegs) {
65948bfdf8aSPeter Brune       ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
66048bfdf8aSPeter Brune     }
66148bfdf8aSPeter Brune   }
66248bfdf8aSPeter Brune   if (fas->downsmooth) {
66348bfdf8aSPeter Brune     ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);
66448bfdf8aSPeter Brune     if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) {
66548bfdf8aSPeter Brune       ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
66648bfdf8aSPeter Brune     }
66748bfdf8aSPeter Brune     if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) {
66848bfdf8aSPeter Brune       ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
66948bfdf8aSPeter Brune     }
67048bfdf8aSPeter Brune    if (snes->ops->computegs && !fas->downsmooth->ops->computegs) {
67148bfdf8aSPeter Brune       ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
67248bfdf8aSPeter Brune     }
6731a266240SBarry Smith   }
67448bfdf8aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
675646217ecSPeter Brune   if (fas->next) {
6766273346dSPeter Brune     if (fas->galerkin) {
6776273346dSPeter Brune       ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr);
6786273346dSPeter Brune     } else {
67948bfdf8aSPeter Brune       if (snes->ops->computefunction && !fas->next->ops->computefunction) {
68048bfdf8aSPeter Brune         ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
68148bfdf8aSPeter Brune       }
68248bfdf8aSPeter Brune       if (snes->ops->computejacobian && !fas->next->ops->computejacobian) {
68348bfdf8aSPeter Brune         ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
68448bfdf8aSPeter Brune       }
68548bfdf8aSPeter Brune       if (snes->ops->computegs && !fas->next->ops->computegs) {
686646217ecSPeter Brune         ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
687646217ecSPeter Brune       }
688646217ecSPeter Brune     }
6896273346dSPeter Brune   }
690421d9b32SPeter Brune   if (snes->dm) {
691421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
692421d9b32SPeter Brune     if (fas->next) {
693421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
694ee78dd50SPeter Brune       if (!fas->next->dm) {
695ee78dd50SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
696ee78dd50SPeter Brune         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
697ee78dd50SPeter Brune       }
698421d9b32SPeter Brune       /* set the interpolation and restriction from the DM */
699ee78dd50SPeter Brune       if (!fas->interpolate) {
700e727c939SJed Brown         ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
701421d9b32SPeter Brune         fas->restrct = fas->interpolate;
702421d9b32SPeter Brune       }
703efe1f98aSPeter Brune       /* set the injection from the DM */
704efe1f98aSPeter Brune       if (!fas->inject) {
705e727c939SJed Brown         ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
706efe1f98aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
707efe1f98aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
708efe1f98aSPeter Brune       }
709421d9b32SPeter Brune     }
710ee78dd50SPeter Brune     /* set the DMs of the pre and post-smoothers here */
711*6bed07a3SBarry Smith     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
712*6bed07a3SBarry Smith     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
713421d9b32SPeter Brune   }
714cc05f883SPeter Brune 
7156273346dSPeter Brune   /* setup FAS work vectors */
7166273346dSPeter Brune   if (fas->galerkin) {
7176273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
7186273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
7196273346dSPeter Brune   }
7206273346dSPeter Brune 
721fa9694d7SPeter Brune   if (fas->next) {
722fa9694d7SPeter Brune    /* gotta set up the solution vector for this to work */
723*6bed07a3SBarry Smith     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
724742fe5e2SPeter Brune     if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);}
725fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
726fa9694d7SPeter Brune   }
727fa9694d7SPeter Brune   /* got to set them all up at once */
728421d9b32SPeter Brune   PetscFunctionReturn(0);
729421d9b32SPeter Brune }
730421d9b32SPeter Brune 
731421d9b32SPeter Brune #undef __FUNCT__
732421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
733421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
734421d9b32SPeter Brune {
735ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
736ee78dd50SPeter Brune   PetscInt       levels = 1;
737d1adcc6fSPeter Brune   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, monflg;
738421d9b32SPeter Brune   PetscErrorCode ierr;
739ee78dd50SPeter Brune   const char     *def_smooth = SNESNRICHARDSON;
740ee78dd50SPeter Brune   char           pre_type[256];
741ee78dd50SPeter Brune   char           post_type[256];
742ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
74307144faaSPeter Brune   SNESFASType    fastype;
744421d9b32SPeter Brune 
745421d9b32SPeter Brune   PetscFunctionBegin;
746c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
747ee78dd50SPeter Brune 
748ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
749ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
750ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
751c732cbdbSBarry Smith     if (!flg && snes->dm) {
752c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
753c732cbdbSBarry Smith       levels++;
754d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
755c732cbdbSBarry Smith     }
756ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
757ee78dd50SPeter Brune   }
758ee78dd50SPeter Brune 
759ee78dd50SPeter Brune   /* type of pre and/or post smoothers -- set both at once */
760ee78dd50SPeter Brune   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
761ee78dd50SPeter Brune   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
76207144faaSPeter Brune   fastype = fas->fastype;
76307144faaSPeter Brune   ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
76407144faaSPeter Brune   if (flg) {
76507144faaSPeter Brune     ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
76607144faaSPeter Brune   }
767d1adcc6fSPeter Brune   ierr = PetscOptionsList("-snes_fas_smoother_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
768d1adcc6fSPeter Brune   if (smoothflg) {
769ee78dd50SPeter Brune     ierr = PetscMemcpy(post_type, pre_type, 256);CHKERRQ(ierr);
770ee78dd50SPeter Brune   } else {
771d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothup_type",  "Nonlinear smoother method","SNESSetType",SNESList,def_smooth,pre_type, 256,&smoothupflg);CHKERRQ(ierr);
772d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_smoothdown_type","Nonlinear smoother method","SNESSetType",SNESList,def_smooth,post_type,256,&smoothdownflg);CHKERRQ(ierr);
773ee78dd50SPeter Brune   }
774ee78dd50SPeter Brune 
775ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
776d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_smoothup","Number of post-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
777d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_smoothdown","Number of pre-smooth iterations","SNESFASSetNumberSmoothUp",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
778d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
779ee78dd50SPeter Brune 
7806273346dSPeter Brune   ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr);
7816273346dSPeter Brune 
782646217ecSPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
783ee78dd50SPeter Brune 
784ee78dd50SPeter Brune   /* other options for the coarsest level */
785ee78dd50SPeter Brune   if (fas->level == 0) {
786d1adcc6fSPeter Brune     ierr = PetscOptionsList("-snes_fas_coarse_smoother_type","Coarsest smoother method","SNESSetType",SNESList,def_smooth,pre_type,256,&smoothflg);CHKERRQ(ierr);
787ee78dd50SPeter Brune   }
788ee78dd50SPeter Brune 
789421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
7908cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
791eff52c0eSPeter Brune 
792d28543b3SPeter Brune   if ((!fas->downsmooth) && ((smoothdownflg || smoothflg) || !snes->usegs)) {
7938cc86e31SPeter Brune     const char     *prefix;
7948cc86e31SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
7958cc86e31SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
7968cc86e31SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
7978cc86e31SPeter Brune     if (fas->level || (fas->levels == 1)) {
798eff52c0eSPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_levels_down_");CHKERRQ(ierr);
7998cc86e31SPeter Brune     } else {
8008cc86e31SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
8018cc86e31SPeter Brune     }
8028cc86e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
8038cc86e31SPeter Brune     ierr = SNESSetType(fas->downsmooth, pre_type);CHKERRQ(ierr);
8048cc86e31SPeter Brune   }
8058cc86e31SPeter Brune 
806d28543b3SPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && ((smoothupflg || smoothflg) || !snes->usegs)) {
80767339d5cSBarry Smith     const char     *prefix;
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);
811eff52c0eSPeter Brune     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_levels_up_");CHKERRQ(ierr);
812293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
813ee78dd50SPeter Brune     ierr = SNESSetType(fas->upsmooth, pre_type);CHKERRQ(ierr);
814ee78dd50SPeter Brune   }
8151a266240SBarry Smith   if (fas->upsmooth) {
8161a266240SBarry Smith     ierr = SNESSetTolerances(fas->upsmooth, 0.0, 0.0, 0.0, fas->max_up_it, 1000);CHKERRQ(ierr);
8171a266240SBarry Smith   }
8181a266240SBarry Smith 
8191a266240SBarry Smith   if (fas->downsmooth) {
8201a266240SBarry Smith     ierr = SNESSetTolerances(fas->downsmooth, 0.0, 0.0, 0.0, fas->max_down_it, 1000);CHKERRQ(ierr);
8211a266240SBarry Smith   }
822ee78dd50SPeter Brune 
823742fe5e2SPeter Brune   if (fas->level != fas->levels - 1) {
824742fe5e2SPeter Brune     ierr = SNESSetTolerances(snes, 0.0, 0.0, 0.0, fas->n_cycles, 1000);CHKERRQ(ierr);
825742fe5e2SPeter Brune   }
826742fe5e2SPeter Brune 
827ee78dd50SPeter Brune   if (monflg) {
828646217ecSPeter Brune     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
829794bee33SPeter Brune     /* set the monitors for the upsmoother and downsmoother */
8302f7ea302SPeter Brune     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
831742fe5e2SPeter Brune     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
832293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
833293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
834d28543b3SPeter Brune   } else {
835d28543b3SPeter Brune     /* unset the monitors on the coarse levels */
836d28543b3SPeter Brune     if (fas->level != fas->levels - 1) {
837d28543b3SPeter Brune       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
838d28543b3SPeter Brune     }
839ee78dd50SPeter Brune   }
840ee78dd50SPeter Brune 
841ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
842d28543b3SPeter Brune   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
843421d9b32SPeter Brune   PetscFunctionReturn(0);
844421d9b32SPeter Brune }
845421d9b32SPeter Brune 
846421d9b32SPeter Brune #undef __FUNCT__
847421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
848421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
849421d9b32SPeter Brune {
850421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
851421d9b32SPeter Brune   PetscBool      iascii;
852421d9b32SPeter Brune   PetscErrorCode ierr;
853421d9b32SPeter Brune   PetscInt levels = fas->levels;
854421d9b32SPeter Brune   PetscInt i;
855421d9b32SPeter Brune 
856421d9b32SPeter Brune   PetscFunctionBegin;
857421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
858421d9b32SPeter Brune   if (iascii) {
859421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
860421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
861421d9b32SPeter Brune     for (i = levels - 1; i >= 0; i--) {
862421d9b32SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
863ee78dd50SPeter Brune       if (fas->upsmooth) {
864ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
865421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
866ee78dd50SPeter Brune         ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
867421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
868421d9b32SPeter Brune       } else {
869ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
870421d9b32SPeter Brune       }
871ee78dd50SPeter Brune       if (fas->downsmooth) {
872ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
873421d9b32SPeter Brune         ierr = PetscViewerASCIIPushTab(viewer);
874ee78dd50SPeter Brune         ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
875421d9b32SPeter Brune         ierr = PetscViewerASCIIPopTab(viewer);
876421d9b32SPeter Brune       } else {
877ee78dd50SPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
878421d9b32SPeter Brune       }
879d28543b3SPeter Brune       if (snes->usegs) {
880eff52c0eSPeter Brune         ierr = PetscViewerASCIIPrintf(viewer, "Using user Gauss-Seidel on level %D\n",  fas->level);CHKERRQ(ierr);
881eff52c0eSPeter Brune       }
882421d9b32SPeter Brune       if (fas->next) fas = (SNES_FAS *)fas->next->data;
883421d9b32SPeter Brune     }
884421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
885421d9b32SPeter Brune   } else {
886421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
887421d9b32SPeter Brune   }
888421d9b32SPeter Brune   PetscFunctionReturn(0);
889421d9b32SPeter Brune }
890421d9b32SPeter Brune 
891421d9b32SPeter Brune #undef __FUNCT__
89239bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth"
89339bd7f45SPeter Brune /*
89439bd7f45SPeter Brune Defines the action of the downsmoother
89539bd7f45SPeter Brune  */
89639bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
89739bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
8982d15c84fSPeter Brune   PetscReal           fnorm;
899742fe5e2SPeter Brune   SNESConvergedReason reason;
90039bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
90139bd7f45SPeter Brune   PetscInt            k;
902421d9b32SPeter Brune   PetscFunctionBegin;
903d1adcc6fSPeter Brune   if (fas->downsmooth) {
904d1adcc6fSPeter Brune     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
905742fe5e2SPeter Brune     /* check convergence reason for the smoother */
906742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
907742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
908742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
909742fe5e2SPeter Brune       PetscFunctionReturn(0);
910742fe5e2SPeter Brune     }
911d28543b3SPeter Brune   } else if (snes->usegs && snes->ops->computegs) {
912794bee33SPeter Brune     if (fas->monitor) {
913794bee33SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
914794bee33SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
915d1adcc6fSPeter Brune       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
916eff52c0eSPeter Brune       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
917d1adcc6fSPeter Brune       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
918794bee33SPeter Brune     }
91907144faaSPeter Brune     for (k = 0; k < fas->max_down_it; k++) {
920646217ecSPeter Brune       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
921cc05f883SPeter Brune       if (fas->monitor) {
922794bee33SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
923794bee33SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
924d1adcc6fSPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
925eff52c0eSPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
926d1adcc6fSPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
927794bee33SPeter Brune       }
928646217ecSPeter Brune     }
929c90fad12SPeter Brune   } else if (snes->pc) {
930c90fad12SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
9312f7ea302SPeter Brune     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
9322f7ea302SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
9332f7ea302SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
9342f7ea302SPeter Brune       PetscFunctionReturn(0);
9352f7ea302SPeter Brune     }
936fe6f9142SPeter Brune   }
93739bd7f45SPeter Brune   PetscFunctionReturn(0);
93839bd7f45SPeter Brune }
93939bd7f45SPeter Brune 
94039bd7f45SPeter Brune 
94139bd7f45SPeter Brune #undef __FUNCT__
94239bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth"
94339bd7f45SPeter Brune /*
94407144faaSPeter Brune Defines the action of the upsmoother
94539bd7f45SPeter Brune  */
94639bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
94739bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
94839bd7f45SPeter Brune   PetscReal           fnorm;
94939bd7f45SPeter Brune   SNESConvergedReason reason;
95039bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
95139bd7f45SPeter Brune   PetscInt            k;
95239bd7f45SPeter Brune   PetscFunctionBegin;
95339bd7f45SPeter Brune   if (fas->upsmooth) {
95439bd7f45SPeter Brune     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
95539bd7f45SPeter Brune     /* check convergence reason for the smoother */
95639bd7f45SPeter Brune     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
95739bd7f45SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
95839bd7f45SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
95939bd7f45SPeter Brune       PetscFunctionReturn(0);
96039bd7f45SPeter Brune     }
96139bd7f45SPeter Brune   } else if (snes->usegs && snes->ops->computegs) {
96239bd7f45SPeter Brune     if (fas->monitor) {
96339bd7f45SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
96439bd7f45SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
96539bd7f45SPeter Brune       ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
96639bd7f45SPeter Brune       ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", 0, fnorm);CHKERRQ(ierr);
96739bd7f45SPeter Brune       ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
96839bd7f45SPeter Brune     }
96907144faaSPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
97039bd7f45SPeter Brune       ierr = SNESComputeGS(snes, B, X);CHKERRQ(ierr);
97139bd7f45SPeter Brune       if (fas->monitor) {
97239bd7f45SPeter Brune         ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
97339bd7f45SPeter Brune         ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
97439bd7f45SPeter Brune         ierr = PetscViewerASCIIAddTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
97539bd7f45SPeter Brune         ierr = PetscViewerASCIIPrintf(fas->monitor, "%d SNES GS Function norm %14.12e\n", k+1, fnorm);CHKERRQ(ierr);
97639bd7f45SPeter Brune         ierr = PetscViewerASCIISubtractTab(fas->monitor,((PetscObject)snes)->tablevel + 2);CHKERRQ(ierr);
97739bd7f45SPeter Brune       }
97839bd7f45SPeter Brune     }
97939bd7f45SPeter Brune   } else if (snes->pc) {
98039bd7f45SPeter Brune     ierr = SNESSolve(snes->pc, B, X);CHKERRQ(ierr);
98139bd7f45SPeter Brune     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
98239bd7f45SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
98339bd7f45SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
98439bd7f45SPeter Brune       PetscFunctionReturn(0);
98539bd7f45SPeter Brune     }
98639bd7f45SPeter Brune   }
98739bd7f45SPeter Brune   PetscFunctionReturn(0);
98839bd7f45SPeter Brune }
98939bd7f45SPeter Brune 
99039bd7f45SPeter Brune #undef __FUNCT__
99139bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection"
99239bd7f45SPeter Brune /*
99339bd7f45SPeter Brune 
99439bd7f45SPeter Brune Performs the FAS coarse correction as:
99539bd7f45SPeter Brune 
99639bd7f45SPeter Brune fine problem: F(x) = 0
99739bd7f45SPeter Brune coarse problem: F^c(x) = b^c
99839bd7f45SPeter Brune 
99939bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
100039bd7f45SPeter Brune 
100139bd7f45SPeter Brune with correction:
100239bd7f45SPeter Brune 
100339bd7f45SPeter Brune 
100439bd7f45SPeter Brune 
100539bd7f45SPeter Brune  */
100639bd7f45SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec B, Vec X, Vec F, Vec X_new) {
100739bd7f45SPeter Brune   PetscErrorCode      ierr;
100839bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
100939bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
101039bd7f45SPeter Brune   SNESConvergedReason reason;
101139bd7f45SPeter Brune   PetscFunctionBegin;
1012fa9694d7SPeter Brune   if (fas->next) {
1013ee78dd50SPeter Brune     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1014794bee33SPeter Brune 
1015c90fad12SPeter Brune     X_c  = fas->next->vec_sol;
1016293a7e31SPeter Brune     Xo_c = fas->next->work[0];
1017c90fad12SPeter Brune     F_c  = fas->next->vec_func;
1018742fe5e2SPeter Brune     B_c  = fas->next->vec_rhs;
1019efe1f98aSPeter Brune 
1020efe1f98aSPeter Brune     /* inject the solution */
1021efe1f98aSPeter Brune     if (fas->inject) {
1022a3cb90a9SPeter Brune       ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
1023efe1f98aSPeter Brune     } else {
1024a3cb90a9SPeter Brune       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
1025a3cb90a9SPeter Brune       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
1026efe1f98aSPeter Brune     }
1027293a7e31SPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1028293a7e31SPeter Brune 
1029293a7e31SPeter Brune     /* restrict the defect */
1030293a7e31SPeter Brune     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1031293a7e31SPeter Brune 
1032c90fad12SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1033ee78dd50SPeter Brune     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1034c90fad12SPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1035742fe5e2SPeter Brune 
1036293a7e31SPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1037c90fad12SPeter Brune 
1038ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
1039ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1040c90fad12SPeter Brune 
1041c90fad12SPeter Brune     /* recurse to the next level */
1042f5a6d4f9SBarry Smith     fas->next->vec_rhs = B_c;
1043742fe5e2SPeter Brune     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
1044742fe5e2SPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1045742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1046742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1047742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1048742fe5e2SPeter Brune       PetscFunctionReturn(0);
1049742fe5e2SPeter Brune     }
1050742fe5e2SPeter Brune     /* fas->next->vec_rhs = PETSC_NULL; */
1051ee78dd50SPeter Brune 
1052fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
1053fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
105439bd7f45SPeter Brune     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr);
1055293a7e31SPeter Brune   }
105639bd7f45SPeter Brune   PetscFunctionReturn(0);
105739bd7f45SPeter Brune }
105839bd7f45SPeter Brune 
105939bd7f45SPeter Brune #undef __FUNCT__
106039bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive"
106139bd7f45SPeter Brune /*
106239bd7f45SPeter Brune 
106339bd7f45SPeter Brune The additive cycle looks like:
106439bd7f45SPeter Brune 
106507144faaSPeter Brune xhat = x
106607144faaSPeter Brune xhat = dS(x, b)
106707144faaSPeter Brune x = coarsecorrection(xhat, b_d)
106807144faaSPeter Brune x = x + nu*(xhat - x);
106939bd7f45SPeter Brune (optional) x = uS(x, b)
107039bd7f45SPeter Brune 
107139bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
107239bd7f45SPeter Brune 
107339bd7f45SPeter Brune  */
107439bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
107507144faaSPeter Brune   Vec                 F, B, Xhat;
107607144faaSPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
107739bd7f45SPeter Brune   PetscErrorCode      ierr;
107807144faaSPeter Brune   SNES_FAS *          fas = (SNES_FAS *)snes->data;
107907144faaSPeter Brune   SNESConvergedReason reason;
108039bd7f45SPeter Brune   PetscFunctionBegin;
108139bd7f45SPeter Brune 
108239bd7f45SPeter Brune   F = snes->vec_func;
108339bd7f45SPeter Brune   B = snes->vec_rhs;
108407144faaSPeter Brune   Xhat = snes->work[1];
108507144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
108607144faaSPeter Brune   /* recurse first */
108707144faaSPeter Brune   if (fas->next) {
108807144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
108939bd7f45SPeter Brune 
109007144faaSPeter Brune     X_c  = fas->next->vec_sol;
109107144faaSPeter Brune     Xo_c = fas->next->work[0];
109207144faaSPeter Brune     F_c  = fas->next->vec_func;
109307144faaSPeter Brune     B_c  = fas->next->vec_rhs;
109439bd7f45SPeter Brune 
109507144faaSPeter Brune     /* inject the solution */
109607144faaSPeter Brune     if (fas->inject) {
109707144faaSPeter Brune       ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr);
109807144faaSPeter Brune     } else {
109907144faaSPeter Brune       ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr);
110007144faaSPeter Brune       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
110107144faaSPeter Brune     }
110207144faaSPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
110307144faaSPeter Brune 
110407144faaSPeter Brune     /* restrict the defect */
110507144faaSPeter Brune     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
110607144faaSPeter Brune 
110707144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
110807144faaSPeter Brune     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
110907144faaSPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
111007144faaSPeter Brune 
111107144faaSPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
111207144faaSPeter Brune 
111307144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
111407144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
111507144faaSPeter Brune 
111607144faaSPeter Brune     /* recurse */
111707144faaSPeter Brune     fas->next->vec_rhs = B_c;
111807144faaSPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
111907144faaSPeter Brune 
112007144faaSPeter Brune     /* smooth on this level */
112107144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
112207144faaSPeter Brune 
112307144faaSPeter Brune    ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
112407144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
112507144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
112607144faaSPeter Brune       PetscFunctionReturn(0);
112707144faaSPeter Brune     }
112807144faaSPeter Brune 
112907144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
113007144faaSPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
113107144faaSPeter Brune     ierr = MatInterpolateAdd(fas->interpolate, X_c, Xhat, Xhat);CHKERRQ(ierr);
113207144faaSPeter Brune 
113307144faaSPeter Brune     /* additive correction */
113407144faaSPeter Brune     ierr = VecAXPY(Xhat, -1.0, X);CHKERRQ(ierr);
113507144faaSPeter Brune     ierr = VecAXPY(X, 0.5, Xhat);CHKERRQ(ierr);
113607144faaSPeter Brune 
113707144faaSPeter Brune   } else {
113807144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
113907144faaSPeter Brune   }
114039bd7f45SPeter Brune   PetscFunctionReturn(0);
114139bd7f45SPeter Brune }
114239bd7f45SPeter Brune 
114339bd7f45SPeter Brune #undef __FUNCT__
114439bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative"
114539bd7f45SPeter Brune /*
114639bd7f45SPeter Brune 
114739bd7f45SPeter Brune Defines the FAS cycle as:
114839bd7f45SPeter Brune 
114939bd7f45SPeter Brune fine problem: F(x) = 0
115039bd7f45SPeter Brune coarse problem: F^c(x) = b^c
115139bd7f45SPeter Brune 
115239bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
115339bd7f45SPeter Brune 
115439bd7f45SPeter Brune correction:
115539bd7f45SPeter Brune 
115639bd7f45SPeter Brune x = x + I(x^c - Rx)
115739bd7f45SPeter Brune 
115839bd7f45SPeter Brune  */
115939bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
116039bd7f45SPeter Brune 
116139bd7f45SPeter Brune   PetscErrorCode      ierr;
116239bd7f45SPeter Brune   Vec                 F,B;
116339bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
116439bd7f45SPeter Brune 
116539bd7f45SPeter Brune   PetscFunctionBegin;
116639bd7f45SPeter Brune   F = snes->vec_func;
116739bd7f45SPeter Brune   B = snes->vec_rhs;
116839bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
116939bd7f45SPeter Brune   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
117039bd7f45SPeter Brune 
117139bd7f45SPeter Brune   ierr = FASCoarseCorrection(snes, B, X, F, X);CHKERRQ(ierr);
117239bd7f45SPeter Brune 
1173c90fad12SPeter Brune   if (fas->level != 0) {
117439bd7f45SPeter Brune     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
1175fe6f9142SPeter Brune   }
1176fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1177fa9694d7SPeter Brune 
1178fa9694d7SPeter Brune   PetscFunctionReturn(0);
1179421d9b32SPeter Brune }
1180421d9b32SPeter Brune 
1181421d9b32SPeter Brune #undef __FUNCT__
1182421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
1183421d9b32SPeter Brune 
1184421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
1185421d9b32SPeter Brune {
1186fa9694d7SPeter Brune   PetscErrorCode ierr;
1187fe6f9142SPeter Brune   PetscInt       i, maxits;
1188f5a6d4f9SBarry Smith   Vec            X, B, F;
1189fe6f9142SPeter Brune   PetscReal      fnorm;
119007144faaSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
1191421d9b32SPeter Brune   PetscFunctionBegin;
1192fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
1193fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1194fa9694d7SPeter Brune   X = snes->vec_sol;
1195fe6f9142SPeter Brune   B = snes->vec_rhs;
1196f5a6d4f9SBarry Smith   F = snes->vec_func;
1197293a7e31SPeter Brune 
1198293a7e31SPeter Brune   /*norm setup */
1199fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1200fe6f9142SPeter Brune   snes->iter = 0;
1201fe6f9142SPeter Brune   snes->norm = 0.;
1202fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1203fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1204fe6f9142SPeter Brune   if (snes->domainerror) {
1205fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1206fe6f9142SPeter Brune     PetscFunctionReturn(0);
1207fe6f9142SPeter Brune   }
1208fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1209fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1210fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1211fe6f9142SPeter Brune   snes->norm = fnorm;
1212fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1213fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1214fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1215fe6f9142SPeter Brune 
1216fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
1217fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
1218fe6f9142SPeter Brune   /* test convergence */
1219fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1220fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
1221fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
1222fe6f9142SPeter Brune     /* Call general purpose update function */
1223646217ecSPeter Brune 
1224fe6f9142SPeter Brune     if (snes->ops->update) {
1225fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1226fe6f9142SPeter Brune     }
122707144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
122839bd7f45SPeter Brune       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
122907144faaSPeter Brune     } else {
123007144faaSPeter Brune       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
123107144faaSPeter Brune     }
1232742fe5e2SPeter Brune 
1233742fe5e2SPeter Brune     /* check for FAS cycle divergence */
1234742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
1235742fe5e2SPeter Brune       PetscFunctionReturn(0);
1236742fe5e2SPeter Brune     }
1237c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1238c90fad12SPeter Brune     /* Monitor convergence */
1239c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1240c90fad12SPeter Brune     snes->iter = i+1;
1241c90fad12SPeter Brune     snes->norm = fnorm;
1242c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1243c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
1244c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1245c90fad12SPeter Brune     /* Test for convergence */
1246c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1247c90fad12SPeter Brune     if (snes->reason) break;
1248fe6f9142SPeter Brune   }
1249fe6f9142SPeter Brune   if (i == maxits) {
1250fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1251fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1252fe6f9142SPeter Brune   }
1253421d9b32SPeter Brune   PetscFunctionReturn(0);
1254421d9b32SPeter Brune }
1255