xref: /petsc/src/snes/impls/fas/fas.c (revision fde0ff241146f0ea9c891c46693ff6f0a90bbe08)
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
7ef536925SPeter Brune 
8ef536925SPeter Brune SNESFAS - Full Approximation Scheme nonlinear multigrid solver.
9421d9b32SPeter Brune 
10421d9b32SPeter Brune The nonlinear problem is solved via the repeated application of nonlinear preconditioners and coarse-grid corrections
11421d9b32SPeter Brune 
12421d9b32SPeter Brune .seealso: SNESCreate(), SNES, SNESSetType(), SNESType (for list of available types)
13421d9b32SPeter Brune M*/
14421d9b32SPeter Brune 
15421d9b32SPeter Brune extern PetscErrorCode SNESDestroy_FAS(SNES snes);
16421d9b32SPeter Brune extern PetscErrorCode SNESSetUp_FAS(SNES snes);
17421d9b32SPeter Brune extern PetscErrorCode SNESSetFromOptions_FAS(SNES snes);
18421d9b32SPeter Brune extern PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer);
19421d9b32SPeter Brune extern PetscErrorCode SNESSolve_FAS(SNES snes);
20421d9b32SPeter Brune extern PetscErrorCode SNESReset_FAS(SNES snes);
216273346dSPeter Brune extern PetscErrorCode SNESFASGalerkinDefaultFunction(SNES, Vec, Vec, void *);
22421d9b32SPeter Brune 
23421d9b32SPeter Brune EXTERN_C_BEGIN
24ddebd997SPeter Brune #undef __FUNCT__
25ddebd997SPeter Brune #define __FUNCT__ "SNESLineSearchSetType_FAS"
26ddebd997SPeter Brune PetscErrorCode  SNESLineSearchSetType_FAS(SNES snes, SNESLineSearchType type)
27ddebd997SPeter Brune {
28ddebd997SPeter Brune   PetscErrorCode ierr;
29ddebd997SPeter Brune   PetscFunctionBegin;
30ddebd997SPeter Brune 
31ddebd997SPeter Brune   switch (type) {
32ddebd997SPeter Brune   case SNES_LS_BASIC:
33ddebd997SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNo,PETSC_NULL);CHKERRQ(ierr);
34ddebd997SPeter Brune     break;
35ddebd997SPeter Brune   case SNES_LS_BASIC_NONORMS:
36ddebd997SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchNoNorms,PETSC_NULL);CHKERRQ(ierr);
37ddebd997SPeter Brune     break;
38ddebd997SPeter Brune   case SNES_LS_QUADRATIC:
39ddebd997SPeter Brune     ierr = SNESLineSearchSet(snes,SNESLineSearchQuadraticSecant,PETSC_NULL);CHKERRQ(ierr);
40ddebd997SPeter Brune     break;
41ddebd997SPeter Brune   default:
42ddebd997SPeter Brune     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP,"Unknown line search type.");
43ddebd997SPeter Brune     break;
44ddebd997SPeter Brune   }
45ddebd997SPeter Brune   snes->ls_type = type;
46ddebd997SPeter Brune   PetscFunctionReturn(0);
47ddebd997SPeter Brune }
48ddebd997SPeter Brune EXTERN_C_END
49ddebd997SPeter Brune 
50ddebd997SPeter Brune EXTERN_C_BEGIN
51421d9b32SPeter Brune 
52421d9b32SPeter Brune #undef __FUNCT__
53421d9b32SPeter Brune #define __FUNCT__ "SNESCreate_FAS"
54421d9b32SPeter Brune PetscErrorCode SNESCreate_FAS(SNES snes)
55421d9b32SPeter Brune {
56421d9b32SPeter Brune   SNES_FAS * fas;
57421d9b32SPeter Brune   PetscErrorCode ierr;
58421d9b32SPeter Brune 
59421d9b32SPeter Brune   PetscFunctionBegin;
60421d9b32SPeter Brune   snes->ops->destroy        = SNESDestroy_FAS;
61421d9b32SPeter Brune   snes->ops->setup          = SNESSetUp_FAS;
62421d9b32SPeter Brune   snes->ops->setfromoptions = SNESSetFromOptions_FAS;
63421d9b32SPeter Brune   snes->ops->view           = SNESView_FAS;
64421d9b32SPeter Brune   snes->ops->solve          = SNESSolve_FAS;
65421d9b32SPeter Brune   snes->ops->reset          = SNESReset_FAS;
66421d9b32SPeter Brune 
67ed020824SBarry Smith   snes->usesksp             = PETSC_FALSE;
68ed020824SBarry Smith   snes->usespc              = PETSC_FALSE;
69ed020824SBarry Smith 
70421d9b32SPeter Brune   ierr = PetscNewLog(snes, SNES_FAS, &fas);CHKERRQ(ierr);
71421d9b32SPeter Brune   snes->data                  = (void*) fas;
72421d9b32SPeter Brune   fas->level                  = 0;
73293a7e31SPeter Brune   fas->levels                 = 1;
74ee78dd50SPeter Brune   fas->n_cycles               = 1;
75ee78dd50SPeter Brune   fas->max_up_it              = 1;
76ee78dd50SPeter Brune   fas->max_down_it            = 1;
77ee78dd50SPeter Brune   fas->upsmooth               = PETSC_NULL;
78ee78dd50SPeter Brune   fas->downsmooth             = PETSC_NULL;
79421d9b32SPeter Brune   fas->next                   = PETSC_NULL;
806273346dSPeter Brune   fas->previous               = PETSC_NULL;
81421d9b32SPeter Brune   fas->interpolate            = PETSC_NULL;
82421d9b32SPeter Brune   fas->restrct                = PETSC_NULL;
83efe1f98aSPeter Brune   fas->inject                 = PETSC_NULL;
84cc05f883SPeter Brune   fas->monitor                = PETSC_NULL;
85cc05f883SPeter Brune   fas->usedmfornumberoflevels = PETSC_FALSE;
86ddebd997SPeter Brune   fas->fastype                = SNES_FAS_MULTIPLICATIVE;
87ddebd997SPeter Brune 
88ddebd997SPeter Brune   ierr = PetscObjectComposeFunctionDynamic((PetscObject)snes,"SNESLineSearchSetType_C","SNESLineSearchSetType_FAS",SNESLineSearchSetType_FAS);CHKERRQ(ierr);
89ddebd997SPeter Brune   ierr = SNESLineSearchSetType(snes, SNES_LS_QUADRATIC);CHKERRQ(ierr);
90ddebd997SPeter Brune 
91421d9b32SPeter Brune   PetscFunctionReturn(0);
92421d9b32SPeter Brune }
93421d9b32SPeter Brune EXTERN_C_END
94421d9b32SPeter Brune 
95421d9b32SPeter Brune #undef __FUNCT__
96421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
97c79ef259SPeter Brune /*@
98722262beSPeter Brune    SNESFASGetLevels - Gets the number of levels in a FAS.
99c79ef259SPeter Brune 
100c79ef259SPeter Brune    Input Parameter:
101c79ef259SPeter Brune .  snes - the preconditioner context
102c79ef259SPeter Brune 
103c79ef259SPeter Brune    Output parameter:
104c79ef259SPeter Brune .  levels - the number of levels
105c79ef259SPeter Brune 
106c79ef259SPeter Brune    Level: advanced
107c79ef259SPeter Brune 
108c79ef259SPeter Brune .keywords: MG, get, levels, multigrid
109c79ef259SPeter Brune 
110c79ef259SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels()
111c79ef259SPeter Brune @*/
112421d9b32SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
113421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
114421d9b32SPeter Brune   PetscFunctionBegin;
115ee1fd11aSPeter Brune   *levels = fas->levels;
116421d9b32SPeter Brune   PetscFunctionReturn(0);
117421d9b32SPeter Brune }
118421d9b32SPeter Brune 
119421d9b32SPeter Brune #undef __FUNCT__
120646217ecSPeter Brune #define __FUNCT__ "SNESFASSetCycles"
121c79ef259SPeter Brune /*@
122c79ef259SPeter Brune    SNESFASSetCycles - Sets the type cycles to use.  Use SNESFASSetCyclesOnLevel() for more
123c79ef259SPeter Brune    complicated cycling.
124c79ef259SPeter Brune 
125c79ef259SPeter Brune    Logically Collective on SNES
126c79ef259SPeter Brune 
127c79ef259SPeter Brune    Input Parameters:
128c79ef259SPeter Brune +  snes   - the multigrid context
129c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
130c79ef259SPeter Brune 
131c79ef259SPeter Brune    Options Database Key:
132c79ef259SPeter Brune $  -snes_fas_cycles 1 or 2
133c79ef259SPeter Brune 
134c79ef259SPeter Brune    Level: advanced
135c79ef259SPeter Brune 
136c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
137c79ef259SPeter Brune 
138c79ef259SPeter Brune .seealso: SNESFASSetCyclesOnLevel()
139c79ef259SPeter Brune @*/
140646217ecSPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
141646217ecSPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
142646217ecSPeter Brune   PetscErrorCode ierr;
143646217ecSPeter Brune   PetscFunctionBegin;
144646217ecSPeter Brune   fas->n_cycles = cycles;
145646217ecSPeter Brune   if (fas->next) {
146646217ecSPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
147646217ecSPeter Brune   }
148646217ecSPeter Brune   PetscFunctionReturn(0);
149646217ecSPeter Brune }
150646217ecSPeter Brune 
151eff52c0eSPeter Brune #undef __FUNCT__
152c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetCyclesOnLevel"
153c79ef259SPeter Brune /*@
154722262beSPeter Brune    SNESFASSetCyclesOnLevel - Sets the type cycles to use on a particular level.
155c79ef259SPeter Brune 
156c79ef259SPeter Brune    Logically Collective on SNES
157c79ef259SPeter Brune 
158c79ef259SPeter Brune    Input Parameters:
159c79ef259SPeter Brune +  snes   - the multigrid context
160c79ef259SPeter Brune .  level  - the level to set the number of cycles on
161c79ef259SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
162c79ef259SPeter Brune 
163c79ef259SPeter Brune    Level: advanced
164c79ef259SPeter Brune 
165c79ef259SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
166c79ef259SPeter Brune 
167c79ef259SPeter Brune .seealso: SNESFASSetCycles()
168c79ef259SPeter Brune @*/
169c79ef259SPeter Brune PetscErrorCode SNESFASSetCyclesOnLevel(SNES snes, PetscInt level, PetscInt cycles) {
170c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
171c79ef259SPeter Brune   PetscInt top_level = fas->level,i;
172c79ef259SPeter Brune 
173c79ef259SPeter Brune   PetscFunctionBegin;
174c79ef259SPeter Brune   if (level > top_level)
175c79ef259SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
176c79ef259SPeter Brune   /* get to the correct level */
177c79ef259SPeter Brune   for (i = fas->level; i > level; i--) {
178c79ef259SPeter Brune     fas = (SNES_FAS *)fas->next->data;
179c79ef259SPeter Brune   }
180c79ef259SPeter Brune   if (fas->level != level)
181c79ef259SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
182c79ef259SPeter Brune   fas->n_cycles = cycles;
183c79ef259SPeter Brune   PetscFunctionReturn(0);
184c79ef259SPeter Brune }
185c79ef259SPeter Brune 
186c79ef259SPeter Brune #undef __FUNCT__
187eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS"
188aeed3662SMatthew G Knepley /*@C
189c79ef259SPeter Brune    SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used.
190c79ef259SPeter Brune    Use SNESFASSetGSOnLevel() for more complicated staging of smoothers
191c79ef259SPeter Brune    and nonlinear preconditioners.
192c79ef259SPeter Brune 
193c79ef259SPeter Brune    Logically Collective on SNES
194c79ef259SPeter Brune 
195c79ef259SPeter Brune    Input Parameters:
196c79ef259SPeter Brune +  snes    - the multigrid context
197c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
198c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
199c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
200c79ef259SPeter Brune 
201c79ef259SPeter Brune    Level: advanced
202c79ef259SPeter Brune 
203c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid
204c79ef259SPeter Brune 
205c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel()
206c79ef259SPeter Brune @*/
207eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
208eff52c0eSPeter Brune   PetscErrorCode ierr = 0;
209eff52c0eSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
210eff52c0eSPeter Brune   PetscFunctionBegin;
211eff52c0eSPeter Brune 
212eff52c0eSPeter Brune   if (gsfunc) {
213eff52c0eSPeter Brune     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
214eff52c0eSPeter Brune     /* push the provided GS up the tree */
215eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
216eff52c0eSPeter Brune   } else if (snes->ops->computegs) {
217eff52c0eSPeter Brune     /* assume that the user has set the GS solver at this level */
218eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
219eff52c0eSPeter Brune   } else if (use_gs) {
220eff52c0eSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %d", fas->level);
221eff52c0eSPeter Brune   }
222eff52c0eSPeter Brune   PetscFunctionReturn(0);
223eff52c0eSPeter Brune }
224eff52c0eSPeter Brune 
225eff52c0eSPeter Brune #undef __FUNCT__
226eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel"
227aeed3662SMatthew G Knepley /*@C
228c79ef259SPeter Brune    SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level.
229c79ef259SPeter Brune 
230c79ef259SPeter Brune    Logically Collective on SNES
231c79ef259SPeter Brune 
232c79ef259SPeter Brune    Input Parameters:
233c79ef259SPeter Brune +  snes    - the multigrid context
234c79ef259SPeter Brune .  level   - the level to set the nonlinear smoother on
235c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
236c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
237c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
238c79ef259SPeter Brune 
239c79ef259SPeter Brune    Level: advanced
240c79ef259SPeter Brune 
241c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
242c79ef259SPeter Brune 
243c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS()
244c79ef259SPeter Brune @*/
245eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
246eff52c0eSPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
247eff52c0eSPeter Brune   PetscErrorCode ierr;
248eff52c0eSPeter Brune   PetscInt       top_level = fas->level,i;
249eff52c0eSPeter Brune   SNES           cur_snes = snes;
250eff52c0eSPeter Brune   PetscFunctionBegin;
251eff52c0eSPeter Brune   if (level > top_level)
252eff52c0eSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetCyclesOnLevel", level);
253eff52c0eSPeter Brune   /* get to the correct level */
254eff52c0eSPeter Brune   for (i = fas->level; i > level; i--) {
255eff52c0eSPeter Brune     fas = (SNES_FAS *)fas->next->data;
256eff52c0eSPeter Brune     cur_snes = fas->next;
257eff52c0eSPeter Brune   }
258eff52c0eSPeter Brune   if (fas->level != level)
259eff52c0eSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetCyclesOnLevel");
260d28543b3SPeter Brune   snes->usegs = use_gs;
261eff52c0eSPeter Brune   if (gsfunc) {
2626273346dSPeter Brune     ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr);
263eff52c0eSPeter Brune   }
264eff52c0eSPeter Brune   PetscFunctionReturn(0);
265eff52c0eSPeter Brune }
266eff52c0eSPeter Brune 
267646217ecSPeter Brune #undef __FUNCT__
268421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
269c79ef259SPeter Brune /*@
270c79ef259SPeter Brune    SNESFASGetSNES - Gets the SNES corresponding to a particular
271c79ef259SPeter Brune    level of the FAS hierarchy.
272c79ef259SPeter Brune 
273c79ef259SPeter Brune    Input Parameters:
274c79ef259SPeter Brune +  snes    - the multigrid context
275c79ef259SPeter Brune    level   - the level to get
276c79ef259SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
277c79ef259SPeter Brune 
278c79ef259SPeter Brune    Level: advanced
279c79ef259SPeter Brune 
280c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
281c79ef259SPeter Brune 
282c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
283c79ef259SPeter Brune @*/
284421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes, PetscInt level, SNES * lsnes) {
285421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
286421d9b32SPeter Brune   PetscInt levels = fas->level;
287421d9b32SPeter Brune   PetscInt i;
288421d9b32SPeter Brune   PetscFunctionBegin;
289421d9b32SPeter Brune   *lsnes = snes;
290421d9b32SPeter Brune   if (fas->level < level) {
291421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES should only be called on a finer SNESFAS instance than the level.");
292421d9b32SPeter Brune   }
293421d9b32SPeter Brune   if (level > levels - 1) {
294421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Level %d doesn't exist in the SNESFAS.");
295421d9b32SPeter Brune   }
296421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
297421d9b32SPeter Brune     *lsnes = fas->next;
298421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
299421d9b32SPeter Brune   }
300421d9b32SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetSNES didn't return the right level!");
301421d9b32SPeter Brune   PetscFunctionReturn(0);
302421d9b32SPeter Brune }
303421d9b32SPeter Brune 
304421d9b32SPeter Brune #undef __FUNCT__
30507144faaSPeter Brune #define __FUNCT__ "SNESFASSetType"
30607144faaSPeter Brune /*@
30707144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS.
30807144faaSPeter Brune e
30907144faaSPeter Brune 
31007144faaSPeter Brune 
31107144faaSPeter Brune @*/
31207144faaSPeter Brune PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
31307144faaSPeter Brune {
31407144faaSPeter Brune   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
31507144faaSPeter Brune 
31607144faaSPeter Brune   PetscFunctionBegin;
31707144faaSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
31807144faaSPeter Brune   PetscValidLogicalCollectiveEnum(snes,fastype,2);
31907144faaSPeter Brune   fas->fastype = fastype;
32007144faaSPeter Brune   PetscFunctionReturn(0);
32107144faaSPeter Brune }
32207144faaSPeter Brune 
32307144faaSPeter Brune 
32407144faaSPeter Brune 
32507144faaSPeter Brune #undef __FUNCT__
326421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
327aeed3662SMatthew G Knepley /*@C
328c79ef259SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
329c79ef259SPeter Brune    Must be called before any other FAS routine.
330c79ef259SPeter Brune 
331c79ef259SPeter Brune    Input Parameters:
332c79ef259SPeter Brune +  snes   - the snes context
333c79ef259SPeter Brune .  levels - the number of levels
334c79ef259SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
335c79ef259SPeter Brune             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
336c79ef259SPeter Brune             Fortran.
337c79ef259SPeter Brune 
338c79ef259SPeter Brune    Level: intermediate
339c79ef259SPeter Brune 
340c79ef259SPeter Brune    Notes:
341c79ef259SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
342c79ef259SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
343c79ef259SPeter Brune 
344c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid
345c79ef259SPeter Brune 
346c79ef259SPeter Brune .seealso: SNESFASGetLevels()
347c79ef259SPeter Brune @*/
348421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
349421d9b32SPeter Brune   PetscErrorCode ierr;
350421d9b32SPeter Brune   PetscInt i;
351421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
3526273346dSPeter Brune   SNES prevsnes;
353421d9b32SPeter Brune   MPI_Comm comm;
354421d9b32SPeter Brune   PetscFunctionBegin;
355ee1fd11aSPeter Brune   comm = ((PetscObject)snes)->comm;
356ee1fd11aSPeter Brune   if (levels == fas->levels) {
357ee1fd11aSPeter Brune     if (!comms)
358ee1fd11aSPeter Brune       PetscFunctionReturn(0);
359ee1fd11aSPeter Brune   }
360421d9b32SPeter Brune   /* user has changed the number of levels; reset */
361421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
362421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
363421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
364ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
3656273346dSPeter Brune   fas->previous = PETSC_NULL;
3666273346dSPeter Brune   prevsnes = snes;
367421d9b32SPeter Brune   /* setup the finest level */
368421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
369421d9b32SPeter Brune     if (comms) comm = comms[i];
370421d9b32SPeter Brune     fas->level = i;
371421d9b32SPeter Brune     fas->levels = levels;
372ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
373421d9b32SPeter Brune     if (i > 0) {
374421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
3754a6b3fb8SBarry Smith       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
376421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
377794bee33SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
3786273346dSPeter Brune       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
3796273346dSPeter Brune       prevsnes = fas->next;
3806273346dSPeter Brune       fas = (SNES_FAS *)prevsnes->data;
381421d9b32SPeter Brune     }
382421d9b32SPeter Brune   }
383421d9b32SPeter Brune   PetscFunctionReturn(0);
384421d9b32SPeter Brune }
385421d9b32SPeter Brune 
386421d9b32SPeter Brune #undef __FUNCT__
387c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp"
388c79ef259SPeter Brune /*@
389c79ef259SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
390c79ef259SPeter Brune    use on all levels.
391c79ef259SPeter Brune 
392fde0ff24SPeter Brune    Logically Collective on SNES
393c79ef259SPeter Brune 
394c79ef259SPeter Brune    Input Parameters:
395c79ef259SPeter Brune +  snes - the multigrid context
396c79ef259SPeter Brune -  n    - the number of smoothing steps
397c79ef259SPeter Brune 
398c79ef259SPeter Brune    Options Database Key:
399d28543b3SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
400c79ef259SPeter Brune 
401c79ef259SPeter Brune    Level: advanced
402c79ef259SPeter Brune 
403c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
404c79ef259SPeter Brune 
405c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
406c79ef259SPeter Brune @*/
407c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
408c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
409c79ef259SPeter Brune   PetscErrorCode ierr = 0;
410c79ef259SPeter Brune   PetscFunctionBegin;
411c79ef259SPeter Brune   fas->max_up_it = n;
412c79ef259SPeter Brune   if (fas->next) {
413c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
414c79ef259SPeter Brune   }
415c79ef259SPeter Brune   PetscFunctionReturn(0);
416c79ef259SPeter Brune }
417c79ef259SPeter Brune 
418c79ef259SPeter Brune #undef __FUNCT__
419c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown"
420c79ef259SPeter Brune /*@
421c79ef259SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
422c79ef259SPeter Brune    use on all levels.
423c79ef259SPeter Brune 
424fde0ff24SPeter Brune    Logically Collective on SNES
425c79ef259SPeter Brune 
426c79ef259SPeter Brune    Input Parameters:
427c79ef259SPeter Brune +  snes - the multigrid context
428c79ef259SPeter Brune -  n    - the number of smoothing steps
429c79ef259SPeter Brune 
430c79ef259SPeter Brune    Options Database Key:
431d28543b3SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
432c79ef259SPeter Brune 
433c79ef259SPeter Brune    Level: advanced
434c79ef259SPeter Brune 
435c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
436c79ef259SPeter Brune 
437c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
438c79ef259SPeter Brune @*/
439c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
440c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
441c79ef259SPeter Brune   PetscErrorCode ierr = 0;
442c79ef259SPeter Brune   PetscFunctionBegin;
443c79ef259SPeter Brune   fas->max_down_it = n;
444c79ef259SPeter Brune   if (fas->next) {
445c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
446c79ef259SPeter Brune   }
447c79ef259SPeter Brune   PetscFunctionReturn(0);
448c79ef259SPeter Brune }
449c79ef259SPeter Brune 
450c79ef259SPeter Brune #undef __FUNCT__
451421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
452c79ef259SPeter Brune /*@
453c79ef259SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
454c79ef259SPeter Brune    interpolation from l-1 to the lth level
455c79ef259SPeter Brune 
456c79ef259SPeter Brune    Input Parameters:
457c79ef259SPeter Brune +  snes      - the multigrid context
458c79ef259SPeter Brune .  mat       - the interpolation operator
459c79ef259SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
460c79ef259SPeter Brune 
461c79ef259SPeter Brune    Level: advanced
462c79ef259SPeter Brune 
463c79ef259SPeter Brune    Notes:
464c79ef259SPeter Brune           Usually this is the same matrix used also to set the restriction
465c79ef259SPeter Brune     for the same level.
466c79ef259SPeter Brune 
467c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
468c79ef259SPeter Brune     out from the matrix size which one it is.
469c79ef259SPeter Brune 
470c79ef259SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
471c79ef259SPeter Brune 
472c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRscale()
473c79ef259SPeter Brune @*/
474421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
475421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
476d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
477421d9b32SPeter Brune 
478421d9b32SPeter Brune   PetscFunctionBegin;
479421d9b32SPeter Brune   if (level > top_level)
480421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInterpolation", level);
481421d9b32SPeter Brune   /* get to the correct level */
482d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
483421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
484421d9b32SPeter Brune   }
485421d9b32SPeter Brune   if (fas->level != level)
486421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInterpolation");
487421d9b32SPeter Brune   fas->interpolate = mat;
488421d9b32SPeter Brune   PetscFunctionReturn(0);
489421d9b32SPeter Brune }
490421d9b32SPeter Brune 
491421d9b32SPeter Brune #undef __FUNCT__
492421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
493c79ef259SPeter Brune /*@
494c79ef259SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
495c79ef259SPeter Brune    from level l to l-1.
496c79ef259SPeter Brune 
497c79ef259SPeter Brune    Input Parameters:
498c79ef259SPeter Brune +  snes  - the multigrid context
499c79ef259SPeter Brune .  mat   - the restriction matrix
500c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
501c79ef259SPeter Brune 
502c79ef259SPeter Brune    Level: advanced
503c79ef259SPeter Brune 
504c79ef259SPeter Brune    Notes:
505c79ef259SPeter Brune           Usually this is the same matrix used also to set the interpolation
506c79ef259SPeter Brune     for the same level.
507c79ef259SPeter Brune 
508c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
509c79ef259SPeter Brune     out from the matrix size which one it is.
510c79ef259SPeter Brune 
511fde0ff24SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
512c79ef259SPeter Brune     is used.
513c79ef259SPeter Brune 
514c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
515c79ef259SPeter Brune 
516c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
517c79ef259SPeter Brune @*/
518421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
519421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
520d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
521421d9b32SPeter Brune 
522421d9b32SPeter Brune   PetscFunctionBegin;
523421d9b32SPeter Brune   if (level > top_level)
524421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
525421d9b32SPeter Brune   /* get to the correct level */
526d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
527421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
528421d9b32SPeter Brune   }
529421d9b32SPeter Brune   if (fas->level != level)
530421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
531421d9b32SPeter Brune   fas->restrct = mat;
532421d9b32SPeter Brune   PetscFunctionReturn(0);
533421d9b32SPeter Brune }
534421d9b32SPeter Brune 
535efe1f98aSPeter Brune 
536421d9b32SPeter Brune #undef __FUNCT__
537421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
538c79ef259SPeter Brune /*@
539c79ef259SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
540c79ef259SPeter Brune    operator from level l to l-1.
541c79ef259SPeter Brune 
542c79ef259SPeter Brune    Input Parameters:
543c79ef259SPeter Brune +  snes   - the multigrid context
544c79ef259SPeter Brune .  rscale - the restriction scaling
545c79ef259SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
546c79ef259SPeter Brune 
547c79ef259SPeter Brune    Level: advanced
548c79ef259SPeter Brune 
549c79ef259SPeter Brune    Notes:
550c79ef259SPeter Brune          This is only used in the case that the injection is not set.
551c79ef259SPeter Brune 
552c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
553c79ef259SPeter Brune 
554c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
555c79ef259SPeter Brune @*/
556421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
557421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
558d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
559421d9b32SPeter Brune 
560421d9b32SPeter Brune   PetscFunctionBegin;
561421d9b32SPeter Brune   if (level > top_level)
562421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetRestriction", level);
563421d9b32SPeter Brune   /* get to the correct level */
564d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
565421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
566421d9b32SPeter Brune   }
567421d9b32SPeter Brune   if (fas->level != level)
568421d9b32SPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetRestriction");
569421d9b32SPeter Brune   fas->rscale = rscale;
570421d9b32SPeter Brune   PetscFunctionReturn(0);
571421d9b32SPeter Brune }
572421d9b32SPeter Brune 
573efe1f98aSPeter Brune 
574efe1f98aSPeter Brune #undef __FUNCT__
575efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection"
576c79ef259SPeter Brune /*@
577c79ef259SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
578c79ef259SPeter Brune    from level l to l-1.
579c79ef259SPeter Brune 
580c79ef259SPeter Brune    Input Parameters:
581c79ef259SPeter Brune +  snes  - the multigrid context
582c79ef259SPeter Brune .  mat   - the restriction matrix
583c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
584c79ef259SPeter Brune 
585c79ef259SPeter Brune    Level: advanced
586c79ef259SPeter Brune 
587c79ef259SPeter Brune    Notes:
588c79ef259SPeter Brune          If you do not set this, the restriction and rscale is used to
589c79ef259SPeter Brune    project the solution instead.
590c79ef259SPeter Brune 
591c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
592c79ef259SPeter Brune 
593c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
594c79ef259SPeter Brune @*/
595efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
596efe1f98aSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
597efe1f98aSPeter Brune   PetscInt top_level = fas->level,i;
598efe1f98aSPeter Brune 
599efe1f98aSPeter Brune   PetscFunctionBegin;
600efe1f98aSPeter Brune   if (level > top_level)
601efe1f98aSPeter Brune     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_OUTOFRANGE, "Bad level number %d in SNESFASSetInjection", level);
602efe1f98aSPeter Brune   /* get to the correct level */
603efe1f98aSPeter Brune   for (i = fas->level; i > level; i--) {
604efe1f98aSPeter Brune     fas = (SNES_FAS *)fas->next->data;
605efe1f98aSPeter Brune   }
606efe1f98aSPeter Brune   if (fas->level != level)
607efe1f98aSPeter Brune     SETERRQ(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "Inconsistent level labelling in SNESFASSetInjection");
608efe1f98aSPeter Brune   fas->inject = mat;
609efe1f98aSPeter Brune   PetscFunctionReturn(0);
610efe1f98aSPeter Brune }
611efe1f98aSPeter Brune 
612421d9b32SPeter Brune #undef __FUNCT__
613421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
614421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
615421d9b32SPeter Brune {
61677df8cc4SPeter Brune   PetscErrorCode ierr = 0;
617421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
618421d9b32SPeter Brune 
619421d9b32SPeter Brune   PetscFunctionBegin;
620742fe5e2SPeter Brune   if (fas->upsmooth)     ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
621742fe5e2SPeter Brune   if (fas->downsmooth)   ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
6223dccd265SPeter Brune   if (fas->inject) {
6233dccd265SPeter Brune     ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
6243dccd265SPeter Brune   }
6253dccd265SPeter Brune   if (fas->interpolate == fas->restrct) {
6263dccd265SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
6273dccd265SPeter Brune     fas->restrct = PETSC_NULL;
6283dccd265SPeter Brune   } else {
6293dccd265SPeter Brune     if (fas->interpolate)  ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
6303dccd265SPeter Brune     if (fas->restrct)      ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
6313dccd265SPeter Brune   }
6323dccd265SPeter Brune   if (fas->rscale)       ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
6333dccd265SPeter Brune 
6343dccd265SPeter Brune   if (fas->upsmooth)   ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
6353dccd265SPeter Brune   if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
636742fe5e2SPeter Brune   if (fas->next)       ierr = SNESReset(fas->next);CHKERRQ(ierr);
637421d9b32SPeter Brune   PetscFunctionReturn(0);
638421d9b32SPeter Brune }
639421d9b32SPeter Brune 
640421d9b32SPeter Brune #undef __FUNCT__
641421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
642421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
643421d9b32SPeter Brune {
644421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
645742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
646421d9b32SPeter Brune 
647421d9b32SPeter Brune   PetscFunctionBegin;
648421d9b32SPeter Brune   /* recursively resets and then destroys */
64951e86f29SPeter Brune   if (fas->upsmooth)     ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
65051e86f29SPeter Brune   if (fas->downsmooth)   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
651421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
652421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
653ee1fd11aSPeter Brune 
654421d9b32SPeter Brune   PetscFunctionReturn(0);
655421d9b32SPeter Brune }
656421d9b32SPeter Brune 
657421d9b32SPeter Brune #undef __FUNCT__
658421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
659421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
660421d9b32SPeter Brune {
66148bfdf8aSPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
662421d9b32SPeter Brune   PetscErrorCode ierr;
663efe1f98aSPeter Brune   VecScatter     injscatter;
664d1adcc6fSPeter Brune   PetscInt       dm_levels;
6653dccd265SPeter Brune   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
666d1adcc6fSPeter Brune 
667421d9b32SPeter Brune   PetscFunctionBegin;
668eff52c0eSPeter Brune 
669cc05f883SPeter Brune   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
670d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
671d1adcc6fSPeter Brune     dm_levels++;
672cc05f883SPeter Brune     if (dm_levels > fas->levels) {
6733dccd265SPeter Brune 
6743dccd265SPeter Brune       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESSetUp_FAS*/
6753dccd265SPeter Brune       vec_sol = snes->vec_sol;
6763dccd265SPeter Brune       vec_func = snes->vec_func;
6773dccd265SPeter Brune       vec_sol_update = snes->vec_sol_update;
6783dccd265SPeter Brune       vec_rhs = snes->vec_rhs;
6793dccd265SPeter Brune       snes->vec_sol = PETSC_NULL;
6803dccd265SPeter Brune       snes->vec_func = PETSC_NULL;
6813dccd265SPeter Brune       snes->vec_sol_update = PETSC_NULL;
6823dccd265SPeter Brune       snes->vec_rhs = PETSC_NULL;
6833dccd265SPeter Brune 
6843dccd265SPeter Brune       /* reset the number of levels */
685d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
686cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
6873dccd265SPeter Brune 
6883dccd265SPeter Brune       snes->vec_sol = vec_sol;
6893dccd265SPeter Brune       snes->vec_func = vec_func;
6903dccd265SPeter Brune       snes->vec_rhs = vec_rhs;
6913dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
692d1adcc6fSPeter Brune     }
693d1adcc6fSPeter Brune   }
694d1adcc6fSPeter Brune 
6953dccd265SPeter Brune   if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
6963dccd265SPeter Brune 
69707144faaSPeter Brune   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
698fde0ff24SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
69907144faaSPeter Brune   } else {
700ddebd997SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
70107144faaSPeter Brune   }
702cc05f883SPeter Brune 
70348bfdf8aSPeter Brune   /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */
70448bfdf8aSPeter Brune   if (fas->upsmooth) {
70548bfdf8aSPeter Brune     if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) {
70648bfdf8aSPeter Brune       ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
70748bfdf8aSPeter Brune     }
70848bfdf8aSPeter Brune     if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) {
70948bfdf8aSPeter Brune       ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
71048bfdf8aSPeter Brune     }
71148bfdf8aSPeter Brune    if (snes->ops->computegs && !fas->upsmooth->ops->computegs) {
71248bfdf8aSPeter Brune       ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
71348bfdf8aSPeter Brune     }
714fde0ff24SPeter Brune    ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);
71548bfdf8aSPeter Brune   }
71648bfdf8aSPeter Brune   if (fas->downsmooth) {
71748bfdf8aSPeter Brune     if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) {
71848bfdf8aSPeter Brune       ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
71948bfdf8aSPeter Brune     }
72048bfdf8aSPeter Brune     if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) {
72148bfdf8aSPeter Brune       ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
72248bfdf8aSPeter Brune     }
72348bfdf8aSPeter Brune    if (snes->ops->computegs && !fas->downsmooth->ops->computegs) {
72448bfdf8aSPeter Brune      ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
72548bfdf8aSPeter Brune     }
726fde0ff24SPeter Brune     ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);
7271a266240SBarry Smith   }
72848bfdf8aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
729646217ecSPeter Brune   if (fas->next) {
7306273346dSPeter Brune     if (fas->galerkin) {
7316273346dSPeter Brune       ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr);
7326273346dSPeter Brune     } else {
73348bfdf8aSPeter Brune       if (snes->ops->computefunction && !fas->next->ops->computefunction) {
73448bfdf8aSPeter Brune         ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
73548bfdf8aSPeter Brune       }
73648bfdf8aSPeter Brune       if (snes->ops->computejacobian && !fas->next->ops->computejacobian) {
73748bfdf8aSPeter Brune         ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
73848bfdf8aSPeter Brune       }
73948bfdf8aSPeter Brune       if (snes->ops->computegs && !fas->next->ops->computegs) {
740646217ecSPeter Brune         ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
741646217ecSPeter Brune       }
742646217ecSPeter Brune     }
7436273346dSPeter Brune   }
744421d9b32SPeter Brune   if (snes->dm) {
745421d9b32SPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
746421d9b32SPeter Brune     if (fas->next) {
747421d9b32SPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
748ee78dd50SPeter Brune       if (!fas->next->dm) {
749ee78dd50SPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
750ee78dd50SPeter Brune         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
751ee78dd50SPeter Brune       }
752421d9b32SPeter Brune       /* set the interpolation and restriction from the DM */
753ee78dd50SPeter Brune       if (!fas->interpolate) {
754e727c939SJed Brown         ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
755421d9b32SPeter Brune         fas->restrct = fas->interpolate;
756421d9b32SPeter Brune       }
757efe1f98aSPeter Brune       /* set the injection from the DM */
758efe1f98aSPeter Brune       if (!fas->inject) {
759e727c939SJed Brown         ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
760efe1f98aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
761efe1f98aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
762efe1f98aSPeter Brune       }
763421d9b32SPeter Brune     }
764ee78dd50SPeter Brune     /* set the DMs of the pre and post-smoothers here */
765*6bed07a3SBarry Smith     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
766*6bed07a3SBarry Smith     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
767421d9b32SPeter Brune   }
768cc05f883SPeter Brune 
7696273346dSPeter Brune   /* setup FAS work vectors */
7706273346dSPeter Brune   if (fas->galerkin) {
7716273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
7726273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
7736273346dSPeter Brune   }
7746273346dSPeter Brune 
775fa9694d7SPeter Brune   if (fas->next) {
776fa9694d7SPeter Brune    /* gotta set up the solution vector for this to work */
777*6bed07a3SBarry Smith     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
778742fe5e2SPeter Brune     if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);}
779fa9694d7SPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
780fa9694d7SPeter Brune   }
781fa9694d7SPeter Brune   /* got to set them all up at once */
782421d9b32SPeter Brune   PetscFunctionReturn(0);
783421d9b32SPeter Brune }
784421d9b32SPeter Brune 
785421d9b32SPeter Brune #undef __FUNCT__
786421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
787421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
788421d9b32SPeter Brune {
789ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
790ee78dd50SPeter Brune   PetscInt       levels = 1;
791fde0ff24SPeter Brune   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg;
792421d9b32SPeter Brune   PetscErrorCode ierr;
793ee78dd50SPeter Brune   const char     *def_smooth = SNESNRICHARDSON;
794ee78dd50SPeter Brune   char           pre_type[256];
795ee78dd50SPeter Brune   char           post_type[256];
796ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
79707144faaSPeter Brune   SNESFASType    fastype;
798fde0ff24SPeter Brune   const char     *optionsprefix;
799fde0ff24SPeter Brune   const char     *prefix;
800421d9b32SPeter Brune 
801421d9b32SPeter Brune   PetscFunctionBegin;
802c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
803ee78dd50SPeter Brune 
804ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
805ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
806ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
807c732cbdbSBarry Smith     if (!flg && snes->dm) {
808c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
809c732cbdbSBarry Smith       levels++;
810d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
811c732cbdbSBarry Smith     }
812ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
813ee78dd50SPeter Brune   }
814ee78dd50SPeter Brune 
815ee78dd50SPeter Brune   /* type of pre and/or post smoothers -- set both at once */
816ee78dd50SPeter Brune   ierr = PetscMemcpy(post_type, def_smooth, 256);CHKERRQ(ierr);
817ee78dd50SPeter Brune   ierr = PetscMemcpy(pre_type, def_smooth, 256);CHKERRQ(ierr);
81807144faaSPeter Brune   fastype = fas->fastype;
81907144faaSPeter Brune   ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
82007144faaSPeter Brune   if (flg) {
82107144faaSPeter Brune     ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
82207144faaSPeter Brune   }
823ee78dd50SPeter Brune 
824fde0ff24SPeter Brune   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
825fde0ff24SPeter Brune 
826fde0ff24SPeter Brune   /* smoother setup options */
827fde0ff24SPeter Brune   ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr);
828fde0ff24SPeter Brune   ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr);
829fde0ff24SPeter Brune   ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr);
830fde0ff24SPeter Brune   if (fas->level == 0) {
831fde0ff24SPeter Brune     ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr);
832fde0ff24SPeter Brune   }
833ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
834fde0ff24SPeter Brune 
835d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
836ee78dd50SPeter Brune 
8376273346dSPeter Brune   ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr);
8386273346dSPeter Brune 
839646217ecSPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
840ee78dd50SPeter Brune 
841ee78dd50SPeter Brune 
842421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
8438cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
844eff52c0eSPeter Brune 
845fde0ff24SPeter Brune   if ((!fas->downsmooth) && smoothcoarseflg) {
846fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
847fde0ff24SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
848fde0ff24SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
849fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
850fde0ff24SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
851fde0ff24SPeter Brune   }
852fde0ff24SPeter Brune 
853fde0ff24SPeter Brune   if ((!fas->downsmooth) && smoothdownflg) {
8548cc86e31SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
8558cc86e31SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
8568cc86e31SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
8578cc86e31SPeter Brune     if (fas->level || (fas->levels == 1)) {
858fde0ff24SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr);
8598cc86e31SPeter Brune     } else {
8608cc86e31SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
8618cc86e31SPeter Brune     }
8628cc86e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
8638cc86e31SPeter Brune   }
8648cc86e31SPeter Brune 
865fde0ff24SPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) {
86667339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
867ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
86867339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
869fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr);
870293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
871ee78dd50SPeter Brune   }
872fde0ff24SPeter Brune 
873fde0ff24SPeter Brune   if ((!fas->downsmooth) && smoothflg) {
874fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
875fde0ff24SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
876fde0ff24SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
877fde0ff24SPeter Brune     if (fas->level || (fas->levels == 1)) {
878fde0ff24SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr);
879fde0ff24SPeter Brune     } else {
880fde0ff24SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
881fde0ff24SPeter Brune     }
882fde0ff24SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
883fde0ff24SPeter Brune   }
884fde0ff24SPeter Brune 
885fde0ff24SPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) {
886fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
887fde0ff24SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
888fde0ff24SPeter Brune     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
889fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr);
890fde0ff24SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
891fde0ff24SPeter Brune   }
892fde0ff24SPeter Brune 
8931a266240SBarry Smith   if (fas->upsmooth) {
894fde0ff24SPeter Brune     ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr);
8951a266240SBarry Smith   }
8961a266240SBarry Smith 
8971a266240SBarry Smith   if (fas->downsmooth) {
898fde0ff24SPeter Brune     ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr);
8991a266240SBarry Smith   }
900ee78dd50SPeter Brune 
901742fe5e2SPeter Brune   if (fas->level != fas->levels - 1) {
902fde0ff24SPeter Brune     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr);
903742fe5e2SPeter Brune   }
904742fe5e2SPeter Brune 
905ee78dd50SPeter Brune   if (monflg) {
906646217ecSPeter Brune     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
907794bee33SPeter Brune     /* set the monitors for the upsmoother and downsmoother */
9082f7ea302SPeter Brune     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
909742fe5e2SPeter Brune     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
910293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
911293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
912d28543b3SPeter Brune   } else {
913d28543b3SPeter Brune     /* unset the monitors on the coarse levels */
914d28543b3SPeter Brune     if (fas->level != fas->levels - 1) {
915d28543b3SPeter Brune       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
916d28543b3SPeter Brune     }
917ee78dd50SPeter Brune   }
918ee78dd50SPeter Brune 
919ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
920d28543b3SPeter Brune   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
921421d9b32SPeter Brune   PetscFunctionReturn(0);
922421d9b32SPeter Brune }
923421d9b32SPeter Brune 
924421d9b32SPeter Brune #undef __FUNCT__
925421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
926421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
927421d9b32SPeter Brune {
928421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
929421d9b32SPeter Brune   PetscBool      iascii;
930421d9b32SPeter Brune   PetscErrorCode ierr;
931421d9b32SPeter Brune 
932421d9b32SPeter Brune   PetscFunctionBegin;
933421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
934421d9b32SPeter Brune   if (iascii) {
935421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %d\n",  fas->levels);CHKERRQ(ierr);
936421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
937421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
938ee78dd50SPeter Brune     if (fas->upsmooth) {
93939a4b097SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
940421d9b32SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);
941ee78dd50SPeter Brune       ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
942421d9b32SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);
943421d9b32SPeter Brune     } else {
944ee78dd50SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
945421d9b32SPeter Brune     }
946ee78dd50SPeter Brune     if (fas->downsmooth) {
94739a4b097SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
948421d9b32SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);
949ee78dd50SPeter Brune       ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
950421d9b32SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);
951421d9b32SPeter Brune     } else {
952ee78dd50SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
953421d9b32SPeter Brune     }
954421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
955421d9b32SPeter Brune   } else {
956421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
957421d9b32SPeter Brune   }
958421d9b32SPeter Brune   PetscFunctionReturn(0);
959421d9b32SPeter Brune }
960421d9b32SPeter Brune 
961421d9b32SPeter Brune #undef __FUNCT__
96239bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth"
96339bd7f45SPeter Brune /*
96439bd7f45SPeter Brune Defines the action of the downsmoother
96539bd7f45SPeter Brune  */
96639bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
96739bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
968fde0ff24SPeter Brune   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
969742fe5e2SPeter Brune   SNESConvergedReason reason;
97039bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
971fde0ff24SPeter Brune   Vec                 G, W, Y;
972fde0ff24SPeter Brune   PetscBool           lssuccess;
973fde0ff24SPeter Brune   PetscInt            k;
974421d9b32SPeter Brune   PetscFunctionBegin;
975fde0ff24SPeter Brune   G = snes->work[1];
976fde0ff24SPeter Brune   W = snes->work[2];
977fde0ff24SPeter Brune   Y = snes->work[3];
978d1adcc6fSPeter Brune   if (fas->downsmooth) {
979d1adcc6fSPeter Brune     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
980742fe5e2SPeter Brune     /* check convergence reason for the smoother */
981742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
982742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
983742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
984742fe5e2SPeter Brune       PetscFunctionReturn(0);
985742fe5e2SPeter Brune     }
986fde0ff24SPeter Brune   } else {
987fde0ff24SPeter Brune     /* basic richardson smoother */
988fde0ff24SPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
989794bee33SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
990794bee33SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
991fde0ff24SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
992fde0ff24SPeter Brune       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
993fde0ff24SPeter Brune       if (!lssuccess) {
994fde0ff24SPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
995fde0ff24SPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
9962f7ea302SPeter Brune           PetscFunctionReturn(0);
9972f7ea302SPeter Brune         }
998fe6f9142SPeter Brune       }
999fde0ff24SPeter Brune       ierr = VecCopy(W, X);CHKERRQ(ierr);
1000fde0ff24SPeter Brune       ierr = VecCopy(G, F);CHKERRQ(ierr);
1001fde0ff24SPeter Brune       fnorm = gnorm;
1002fde0ff24SPeter Brune     }
1003fde0ff24SPeter Brune   }
100439bd7f45SPeter Brune   PetscFunctionReturn(0);
100539bd7f45SPeter Brune }
100639bd7f45SPeter Brune 
100739bd7f45SPeter Brune 
100839bd7f45SPeter Brune #undef __FUNCT__
100939bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth"
101039bd7f45SPeter Brune /*
101107144faaSPeter Brune Defines the action of the upsmoother
101239bd7f45SPeter Brune  */
101339bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
101439bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
1015fde0ff24SPeter Brune   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
101639bd7f45SPeter Brune   SNESConvergedReason reason;
101739bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
1018fde0ff24SPeter Brune   Vec                 G, W, Y;
1019fde0ff24SPeter Brune   PetscBool           lssuccess;
1020fde0ff24SPeter Brune   PetscInt            k;
102139bd7f45SPeter Brune   PetscFunctionBegin;
1022fde0ff24SPeter Brune   G = snes->work[1];
1023fde0ff24SPeter Brune   W = snes->work[2];
1024fde0ff24SPeter Brune   Y = snes->work[3];
102539bd7f45SPeter Brune   if (fas->upsmooth) {
1026fde0ff24SPeter Brune     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
102739bd7f45SPeter Brune     /* check convergence reason for the smoother */
1028fde0ff24SPeter Brune     ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr);
102939bd7f45SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
103039bd7f45SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
103139bd7f45SPeter Brune       PetscFunctionReturn(0);
103239bd7f45SPeter Brune     }
1033fde0ff24SPeter Brune   } else {
1034fde0ff24SPeter Brune     /* basic richardson smoother */
1035fde0ff24SPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
103639bd7f45SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
103739bd7f45SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1038fde0ff24SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
1039fde0ff24SPeter Brune       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
1040fde0ff24SPeter Brune       if (!lssuccess) {
1041fde0ff24SPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
1042fde0ff24SPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
104339bd7f45SPeter Brune           PetscFunctionReturn(0);
104439bd7f45SPeter Brune         }
104539bd7f45SPeter Brune       }
1046fde0ff24SPeter Brune       ierr = VecCopy(W, X);CHKERRQ(ierr);
1047fde0ff24SPeter Brune       ierr = VecCopy(G, F);CHKERRQ(ierr);
1048fde0ff24SPeter Brune       fnorm = gnorm;
1049fde0ff24SPeter Brune     }
1050fde0ff24SPeter Brune   }
105139bd7f45SPeter Brune   PetscFunctionReturn(0);
105239bd7f45SPeter Brune }
105339bd7f45SPeter Brune 
105439bd7f45SPeter Brune #undef __FUNCT__
105539bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection"
105639bd7f45SPeter Brune /*
105739bd7f45SPeter Brune 
105839bd7f45SPeter Brune Performs the FAS coarse correction as:
105939bd7f45SPeter Brune 
106039bd7f45SPeter Brune fine problem: F(x) = 0
106139bd7f45SPeter Brune coarse problem: F^c(x) = b^c
106239bd7f45SPeter Brune 
106339bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
106439bd7f45SPeter Brune 
106539bd7f45SPeter Brune with correction:
106639bd7f45SPeter Brune 
106739bd7f45SPeter Brune 
106839bd7f45SPeter Brune 
106939bd7f45SPeter Brune  */
107039a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
107139bd7f45SPeter Brune   PetscErrorCode      ierr;
107239bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
107339bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
107439bd7f45SPeter Brune   SNESConvergedReason reason;
107539bd7f45SPeter Brune   PetscFunctionBegin;
1076fa9694d7SPeter Brune   if (fas->next) {
1077ee78dd50SPeter Brune     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1078794bee33SPeter Brune 
1079c90fad12SPeter Brune     X_c  = fas->next->vec_sol;
1080293a7e31SPeter Brune     Xo_c = fas->next->work[0];
1081c90fad12SPeter Brune     F_c  = fas->next->vec_func;
1082742fe5e2SPeter Brune     B_c  = fas->next->vec_rhs;
1083efe1f98aSPeter Brune 
1084efe1f98aSPeter Brune     /* inject the solution */
1085efe1f98aSPeter Brune     if (fas->inject) {
1086a3cb90a9SPeter Brune       ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
1087efe1f98aSPeter Brune     } else {
1088a3cb90a9SPeter Brune       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
1089a3cb90a9SPeter Brune       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
1090efe1f98aSPeter Brune     }
1091293a7e31SPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1092293a7e31SPeter Brune 
1093293a7e31SPeter Brune     /* restrict the defect */
1094293a7e31SPeter Brune     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1095293a7e31SPeter Brune 
1096c90fad12SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1097ee78dd50SPeter Brune     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1098c90fad12SPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1099742fe5e2SPeter Brune 
1100293a7e31SPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1101c90fad12SPeter Brune 
1102ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
1103ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1104c90fad12SPeter Brune 
1105c90fad12SPeter Brune     /* recurse to the next level */
1106f5a6d4f9SBarry Smith     fas->next->vec_rhs = B_c;
1107742fe5e2SPeter Brune     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
1108742fe5e2SPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1109742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1110742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1111742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1112742fe5e2SPeter Brune       PetscFunctionReturn(0);
1113742fe5e2SPeter Brune     }
1114742fe5e2SPeter Brune     /* fas->next->vec_rhs = PETSC_NULL; */
1115ee78dd50SPeter Brune 
1116fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
1117fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
111839bd7f45SPeter Brune     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr);
1119293a7e31SPeter Brune   }
112039bd7f45SPeter Brune   PetscFunctionReturn(0);
112139bd7f45SPeter Brune }
112239bd7f45SPeter Brune 
112339bd7f45SPeter Brune #undef __FUNCT__
112439bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive"
112539bd7f45SPeter Brune /*
112639bd7f45SPeter Brune 
112739bd7f45SPeter Brune The additive cycle looks like:
112839bd7f45SPeter Brune 
112907144faaSPeter Brune xhat = x
113007144faaSPeter Brune xhat = dS(x, b)
113107144faaSPeter Brune x = coarsecorrection(xhat, b_d)
113207144faaSPeter Brune x = x + nu*(xhat - x);
113339bd7f45SPeter Brune (optional) x = uS(x, b)
113439bd7f45SPeter Brune 
113539bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
113639bd7f45SPeter Brune 
113739bd7f45SPeter Brune  */
113839bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
113907144faaSPeter Brune   Vec                 F, B, Xhat;
1140ddebd997SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c, G, W;
114139bd7f45SPeter Brune   PetscErrorCode      ierr;
114207144faaSPeter Brune   SNES_FAS *          fas = (SNES_FAS *)snes->data;
114307144faaSPeter Brune   SNESConvergedReason reason;
1144ddebd997SPeter Brune   PetscReal           xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.;
1145ddebd997SPeter Brune   PetscBool           lssucceed;
114639bd7f45SPeter Brune   PetscFunctionBegin;
114739bd7f45SPeter Brune 
114839bd7f45SPeter Brune   F = snes->vec_func;
114939bd7f45SPeter Brune   B = snes->vec_rhs;
1150fde0ff24SPeter Brune   Xhat = snes->work[3];
1151fde0ff24SPeter Brune   G    = snes->work[1];
1152fde0ff24SPeter Brune   W    = snes->work[2];
115307144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
115407144faaSPeter Brune   /* recurse first */
115507144faaSPeter Brune   if (fas->next) {
115607144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
115739bd7f45SPeter Brune 
115807144faaSPeter Brune     X_c  = fas->next->vec_sol;
115907144faaSPeter Brune     Xo_c = fas->next->work[0];
116007144faaSPeter Brune     F_c  = fas->next->vec_func;
116107144faaSPeter Brune     B_c  = fas->next->vec_rhs;
116239bd7f45SPeter Brune 
116307144faaSPeter Brune     /* inject the solution */
116407144faaSPeter Brune     if (fas->inject) {
116507144faaSPeter Brune       ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr);
116607144faaSPeter Brune     } else {
116707144faaSPeter Brune       ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr);
116807144faaSPeter Brune       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
116907144faaSPeter Brune     }
117007144faaSPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
117107144faaSPeter Brune 
117207144faaSPeter Brune     /* restrict the defect */
117307144faaSPeter Brune     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
117407144faaSPeter Brune 
117507144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
117607144faaSPeter Brune     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
117707144faaSPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
117807144faaSPeter Brune 
117907144faaSPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
118007144faaSPeter Brune 
118107144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
118207144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
118307144faaSPeter Brune 
118407144faaSPeter Brune     /* recurse */
118507144faaSPeter Brune     fas->next->vec_rhs = B_c;
118607144faaSPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
118707144faaSPeter Brune 
118807144faaSPeter Brune     /* smooth on this level */
118907144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
119007144faaSPeter Brune 
119107144faaSPeter Brune    ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
119207144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
119307144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
119407144faaSPeter Brune       PetscFunctionReturn(0);
119507144faaSPeter Brune     }
119607144faaSPeter Brune 
119707144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
119807144faaSPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1199ddebd997SPeter Brune     ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr);
120007144faaSPeter Brune 
1201ddebd997SPeter Brune     /* additive correction of the coarse direction*/
1202ddebd997SPeter Brune     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1203ddebd997SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1204eb1825c3SPeter Brune     ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr);
1205ddebd997SPeter Brune     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1206ddebd997SPeter Brune     ierr = VecCopy(W, X);CHKERRQ(ierr);
1207ddebd997SPeter Brune     ierr = VecCopy(G, F);CHKERRQ(ierr);
1208ddebd997SPeter Brune     fnorm = gnorm;
120907144faaSPeter Brune   } else {
121007144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
121107144faaSPeter Brune   }
121239bd7f45SPeter Brune   PetscFunctionReturn(0);
121339bd7f45SPeter Brune }
121439bd7f45SPeter Brune 
121539bd7f45SPeter Brune #undef __FUNCT__
121639bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative"
121739bd7f45SPeter Brune /*
121839bd7f45SPeter Brune 
121939bd7f45SPeter Brune Defines the FAS cycle as:
122039bd7f45SPeter Brune 
122139bd7f45SPeter Brune fine problem: F(x) = 0
122239bd7f45SPeter Brune coarse problem: F^c(x) = b^c
122339bd7f45SPeter Brune 
122439bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
122539bd7f45SPeter Brune 
122639bd7f45SPeter Brune correction:
122739bd7f45SPeter Brune 
122839bd7f45SPeter Brune x = x + I(x^c - Rx)
122939bd7f45SPeter Brune 
123039bd7f45SPeter Brune  */
123139bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
123239bd7f45SPeter Brune 
123339bd7f45SPeter Brune   PetscErrorCode      ierr;
123439bd7f45SPeter Brune   Vec                 F,B;
123539bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
123639bd7f45SPeter Brune 
123739bd7f45SPeter Brune   PetscFunctionBegin;
123839bd7f45SPeter Brune   F = snes->vec_func;
123939bd7f45SPeter Brune   B = snes->vec_rhs;
124039bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
124139bd7f45SPeter Brune   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
124239bd7f45SPeter Brune 
124339a4b097SPeter Brune   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
124439bd7f45SPeter Brune 
1245c90fad12SPeter Brune   if (fas->level != 0) {
124639bd7f45SPeter Brune     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
1247fe6f9142SPeter Brune   }
1248fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1249fa9694d7SPeter Brune 
1250fa9694d7SPeter Brune   PetscFunctionReturn(0);
1251421d9b32SPeter Brune }
1252421d9b32SPeter Brune 
1253421d9b32SPeter Brune #undef __FUNCT__
1254421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
1255421d9b32SPeter Brune 
1256421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
1257421d9b32SPeter Brune {
1258fa9694d7SPeter Brune   PetscErrorCode ierr;
1259fe6f9142SPeter Brune   PetscInt       i, maxits;
1260ddb5aff1SPeter Brune   Vec            X, F;
1261fe6f9142SPeter Brune   PetscReal      fnorm;
126207144faaSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
1263421d9b32SPeter Brune   PetscFunctionBegin;
1264fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
1265fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1266fa9694d7SPeter Brune   X = snes->vec_sol;
1267f5a6d4f9SBarry Smith   F = snes->vec_func;
1268293a7e31SPeter Brune 
1269293a7e31SPeter Brune   /*norm setup */
1270fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1271fe6f9142SPeter Brune   snes->iter = 0;
1272fe6f9142SPeter Brune   snes->norm = 0.;
1273fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1274fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1275fe6f9142SPeter Brune   if (snes->domainerror) {
1276fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1277fe6f9142SPeter Brune     PetscFunctionReturn(0);
1278fe6f9142SPeter Brune   }
1279fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1280fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1281fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1282fe6f9142SPeter Brune   snes->norm = fnorm;
1283fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1284fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1285fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1286fe6f9142SPeter Brune 
1287fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
1288fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
1289fe6f9142SPeter Brune   /* test convergence */
1290fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1291fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
1292fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
1293fe6f9142SPeter Brune     /* Call general purpose update function */
1294646217ecSPeter Brune 
1295fe6f9142SPeter Brune     if (snes->ops->update) {
1296fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1297fe6f9142SPeter Brune     }
129807144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
129939bd7f45SPeter Brune       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
130007144faaSPeter Brune     } else {
130107144faaSPeter Brune       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
130207144faaSPeter Brune     }
1303742fe5e2SPeter Brune 
1304742fe5e2SPeter Brune     /* check for FAS cycle divergence */
1305742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
1306742fe5e2SPeter Brune       PetscFunctionReturn(0);
1307742fe5e2SPeter Brune     }
1308c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1309c90fad12SPeter Brune     /* Monitor convergence */
1310c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1311c90fad12SPeter Brune     snes->iter = i+1;
1312c90fad12SPeter Brune     snes->norm = fnorm;
1313c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1314c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
1315c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1316c90fad12SPeter Brune     /* Test for convergence */
1317c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1318c90fad12SPeter Brune     if (snes->reason) break;
1319fe6f9142SPeter Brune   }
1320fe6f9142SPeter Brune   if (i == maxits) {
1321fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1322fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1323fe6f9142SPeter Brune   }
1324421d9b32SPeter Brune   PetscFunctionReturn(0);
1325421d9b32SPeter Brune }
1326