xref: /petsc/src/snes/impls/fas/fas.c (revision 2e8ce2484af3a80acdce8ced464caed9dd1fa784)
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 /*@
98*2e8ce248SJed Brown    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
99c79ef259SPeter Brune 
100c79ef259SPeter Brune    Input Parameter:
101*2e8ce248SJed Brown .  snes - the nonlinear solver 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 /*@
122*2e8ce248SJed Brown    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  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;
174*2e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D cycle count from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
175c79ef259SPeter Brune   /* get to the correct level */
176c79ef259SPeter Brune   for (i = fas->level; i > level; i--) {
177c79ef259SPeter Brune     fas = (SNES_FAS *)fas->next->data;
178c79ef259SPeter Brune   }
179*2e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
180c79ef259SPeter Brune   fas->n_cycles = cycles;
181c79ef259SPeter Brune   PetscFunctionReturn(0);
182c79ef259SPeter Brune }
183c79ef259SPeter Brune 
184c79ef259SPeter Brune #undef __FUNCT__
185eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGS"
186aeed3662SMatthew G Knepley /*@C
187c79ef259SPeter Brune    SNESFASSetGS - Sets a nonlinear GS smoother and if it should be used.
188c79ef259SPeter Brune    Use SNESFASSetGSOnLevel() for more complicated staging of smoothers
189c79ef259SPeter Brune    and nonlinear preconditioners.
190c79ef259SPeter Brune 
191c79ef259SPeter Brune    Logically Collective on SNES
192c79ef259SPeter Brune 
193c79ef259SPeter Brune    Input Parameters:
194c79ef259SPeter Brune +  snes    - the multigrid context
195c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
196c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
197c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
198c79ef259SPeter Brune 
199c79ef259SPeter Brune    Level: advanced
200c79ef259SPeter Brune 
201c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, gauss-seidel, multigrid
202c79ef259SPeter Brune 
203c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGSOnLevel()
204c79ef259SPeter Brune @*/
205eff52c0eSPeter Brune PetscErrorCode SNESFASSetGS(SNES snes, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
206eff52c0eSPeter Brune   PetscErrorCode ierr = 0;
207eff52c0eSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
208eff52c0eSPeter Brune   PetscFunctionBegin;
209eff52c0eSPeter Brune 
210eff52c0eSPeter Brune   if (gsfunc) {
211eff52c0eSPeter Brune     ierr = SNESSetGS(snes, gsfunc, ctx);CHKERRQ(ierr);
212eff52c0eSPeter Brune     /* push the provided GS up the tree */
213eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, gsfunc, ctx, use_gs);CHKERRQ(ierr);
214eff52c0eSPeter Brune   } else if (snes->ops->computegs) {
215eff52c0eSPeter Brune     /* assume that the user has set the GS solver at this level */
216eff52c0eSPeter Brune     if (fas->next) ierr = SNESFASSetGS(fas->next, PETSC_NULL, PETSC_NULL, use_gs);CHKERRQ(ierr);
217eff52c0eSPeter Brune   } else if (use_gs) {
218*2e8ce248SJed Brown     SETERRQ1(((PetscObject)snes)->comm, PETSC_ERR_ARG_WRONG, "No user Gauss-Seidel function provided in SNESFASSetGS on level %D", fas->level);
219eff52c0eSPeter Brune   }
220eff52c0eSPeter Brune   PetscFunctionReturn(0);
221eff52c0eSPeter Brune }
222eff52c0eSPeter Brune 
223eff52c0eSPeter Brune #undef __FUNCT__
224eff52c0eSPeter Brune #define __FUNCT__ "SNESFASSetGSOnLevel"
225aeed3662SMatthew G Knepley /*@C
226c79ef259SPeter Brune    SNESFASSetGSOnLevel - Sets the nonlinear smoother on a particular level.
227c79ef259SPeter Brune 
228c79ef259SPeter Brune    Logically Collective on SNES
229c79ef259SPeter Brune 
230c79ef259SPeter Brune    Input Parameters:
231c79ef259SPeter Brune +  snes    - the multigrid context
232c79ef259SPeter Brune .  level   - the level to set the nonlinear smoother on
233c79ef259SPeter Brune .  gsfunc  - the nonlinear smoother function
234c79ef259SPeter Brune .  ctx     - the user context for the nonlinear smoother
235c79ef259SPeter Brune -  use_gs  - whether to use the nonlinear smoother or not
236c79ef259SPeter Brune 
237c79ef259SPeter Brune    Level: advanced
238c79ef259SPeter Brune 
239c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
240c79ef259SPeter Brune 
241c79ef259SPeter Brune .seealso: SNESSetGS(), SNESFASSetGS()
242c79ef259SPeter Brune @*/
243eff52c0eSPeter Brune PetscErrorCode SNESFASSetGSOnLevel(SNES snes, PetscInt level, PetscErrorCode (*gsfunc)(SNES,Vec,Vec,void *), void * ctx, PetscBool use_gs) {
244eff52c0eSPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
245eff52c0eSPeter Brune   PetscErrorCode ierr;
246eff52c0eSPeter Brune   PetscInt       top_level = fas->level,i;
247eff52c0eSPeter Brune   SNES           cur_snes = snes;
248eff52c0eSPeter Brune   PetscFunctionBegin;
249*2e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D GS from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
250eff52c0eSPeter Brune   /* get to the correct level */
251eff52c0eSPeter Brune   for (i = fas->level; i > level; i--) {
252eff52c0eSPeter Brune     fas = (SNES_FAS *)fas->next->data;
253eff52c0eSPeter Brune     cur_snes = fas->next;
254eff52c0eSPeter Brune   }
255*2e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
256eff52c0eSPeter Brune   if (gsfunc) {
2576273346dSPeter Brune     ierr = SNESSetGS(cur_snes, gsfunc, ctx);CHKERRQ(ierr);
258eff52c0eSPeter Brune   }
259eff52c0eSPeter Brune   PetscFunctionReturn(0);
260eff52c0eSPeter Brune }
261eff52c0eSPeter Brune 
262646217ecSPeter Brune #undef __FUNCT__
263421d9b32SPeter Brune #define __FUNCT__ "SNESFASGetSNES"
264c79ef259SPeter Brune /*@
265c79ef259SPeter Brune    SNESFASGetSNES - Gets the SNES corresponding to a particular
266c79ef259SPeter Brune    level of the FAS hierarchy.
267c79ef259SPeter Brune 
268c79ef259SPeter Brune    Input Parameters:
269c79ef259SPeter Brune +  snes    - the multigrid context
270c79ef259SPeter Brune    level   - the level to get
271c79ef259SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
272c79ef259SPeter Brune 
273c79ef259SPeter Brune    Level: advanced
274c79ef259SPeter Brune 
275c79ef259SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
276c79ef259SPeter Brune 
277c79ef259SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
278c79ef259SPeter Brune @*/
279421d9b32SPeter Brune PetscErrorCode SNESFASGetSNES(SNES snes,PetscInt level,SNES *lsnes) {
280421d9b32SPeter Brune   SNES_FAS *fas = (SNES_FAS*)snes->data;
281421d9b32SPeter Brune   PetscInt i;
282*2e8ce248SJed Brown 
283421d9b32SPeter Brune   PetscFunctionBegin;
284*2e8ce248SJed Brown   if (level > fas->levels-1) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
285*2e8ce248SJed Brown   if (fas->level < level) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS on level %D, must call from SNES at least as fine as level",level,fas->level);
286421d9b32SPeter Brune   *lsnes = snes;
287421d9b32SPeter Brune   for (i = fas->level; i > level; i--) {
288421d9b32SPeter Brune     *lsnes = fas->next;
289421d9b32SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
290421d9b32SPeter Brune   }
291*2e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
292421d9b32SPeter Brune   PetscFunctionReturn(0);
293421d9b32SPeter Brune }
294421d9b32SPeter Brune 
295421d9b32SPeter Brune #undef __FUNCT__
29607144faaSPeter Brune #define __FUNCT__ "SNESFASSetType"
29707144faaSPeter Brune /*@
29807144faaSPeter Brune SNESFASSetType - Sets the update and correction type used for FAS.
29907144faaSPeter Brune e
30007144faaSPeter Brune 
30107144faaSPeter Brune 
30207144faaSPeter Brune @*/
30307144faaSPeter Brune PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
30407144faaSPeter Brune {
30507144faaSPeter Brune   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
30607144faaSPeter Brune 
30707144faaSPeter Brune   PetscFunctionBegin;
30807144faaSPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
30907144faaSPeter Brune   PetscValidLogicalCollectiveEnum(snes,fastype,2);
31007144faaSPeter Brune   fas->fastype = fastype;
31107144faaSPeter Brune   PetscFunctionReturn(0);
31207144faaSPeter Brune }
31307144faaSPeter Brune 
31407144faaSPeter Brune #undef __FUNCT__
315421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
316aeed3662SMatthew G Knepley /*@C
317c79ef259SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
318c79ef259SPeter Brune    Must be called before any other FAS routine.
319c79ef259SPeter Brune 
320c79ef259SPeter Brune    Input Parameters:
321c79ef259SPeter Brune +  snes   - the snes context
322c79ef259SPeter Brune .  levels - the number of levels
323c79ef259SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
324c79ef259SPeter Brune             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
325c79ef259SPeter Brune             Fortran.
326c79ef259SPeter Brune 
327c79ef259SPeter Brune    Level: intermediate
328c79ef259SPeter Brune 
329c79ef259SPeter Brune    Notes:
330c79ef259SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
331c79ef259SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
332c79ef259SPeter Brune 
333c79ef259SPeter Brune .keywords: FAS, MG, set, levels, multigrid
334c79ef259SPeter Brune 
335c79ef259SPeter Brune .seealso: SNESFASGetLevels()
336c79ef259SPeter Brune @*/
337421d9b32SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
338421d9b32SPeter Brune   PetscErrorCode ierr;
339421d9b32SPeter Brune   PetscInt i;
340421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
3416273346dSPeter Brune   SNES prevsnes;
342421d9b32SPeter Brune   MPI_Comm comm;
343421d9b32SPeter Brune   PetscFunctionBegin;
344ee1fd11aSPeter Brune   comm = ((PetscObject)snes)->comm;
345ee1fd11aSPeter Brune   if (levels == fas->levels) {
346ee1fd11aSPeter Brune     if (!comms)
347ee1fd11aSPeter Brune       PetscFunctionReturn(0);
348ee1fd11aSPeter Brune   }
349421d9b32SPeter Brune   /* user has changed the number of levels; reset */
350421d9b32SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
351421d9b32SPeter Brune   /* destroy any coarser levels if necessary */
352421d9b32SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
353ee1fd11aSPeter Brune   fas->next = PETSC_NULL;
3546273346dSPeter Brune   fas->previous = PETSC_NULL;
3556273346dSPeter Brune   prevsnes = snes;
356421d9b32SPeter Brune   /* setup the finest level */
357421d9b32SPeter Brune   for (i = levels - 1; i >= 0; i--) {
358421d9b32SPeter Brune     if (comms) comm = comms[i];
359421d9b32SPeter Brune     fas->level = i;
360421d9b32SPeter Brune     fas->levels = levels;
361ee1fd11aSPeter Brune     fas->next = PETSC_NULL;
362421d9b32SPeter Brune     if (i > 0) {
363421d9b32SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
3644a6b3fb8SBarry Smith       ierr = SNESSetOptionsPrefix(fas->next,((PetscObject)snes)->prefix);CHKERRQ(ierr);
365421d9b32SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
366794bee33SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
3676273346dSPeter Brune       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
3686273346dSPeter Brune       prevsnes = fas->next;
3696273346dSPeter Brune       fas = (SNES_FAS *)prevsnes->data;
370421d9b32SPeter Brune     }
371421d9b32SPeter Brune   }
372421d9b32SPeter Brune   PetscFunctionReturn(0);
373421d9b32SPeter Brune }
374421d9b32SPeter Brune 
375421d9b32SPeter Brune #undef __FUNCT__
376c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp"
377c79ef259SPeter Brune /*@
378c79ef259SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
379c79ef259SPeter Brune    use on all levels.
380c79ef259SPeter Brune 
381fde0ff24SPeter Brune    Logically Collective on SNES
382c79ef259SPeter Brune 
383c79ef259SPeter Brune    Input Parameters:
384c79ef259SPeter Brune +  snes - the multigrid context
385c79ef259SPeter Brune -  n    - the number of smoothing steps
386c79ef259SPeter Brune 
387c79ef259SPeter Brune    Options Database Key:
388d28543b3SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
389c79ef259SPeter Brune 
390c79ef259SPeter Brune    Level: advanced
391c79ef259SPeter Brune 
392c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
393c79ef259SPeter Brune 
394c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
395c79ef259SPeter Brune @*/
396c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
397c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
398c79ef259SPeter Brune   PetscErrorCode ierr = 0;
399c79ef259SPeter Brune   PetscFunctionBegin;
400c79ef259SPeter Brune   fas->max_up_it = n;
401c79ef259SPeter Brune   if (fas->next) {
402c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
403c79ef259SPeter Brune   }
404c79ef259SPeter Brune   PetscFunctionReturn(0);
405c79ef259SPeter Brune }
406c79ef259SPeter Brune 
407c79ef259SPeter Brune #undef __FUNCT__
408c79ef259SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown"
409c79ef259SPeter Brune /*@
410c79ef259SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
411c79ef259SPeter Brune    use on all levels.
412c79ef259SPeter Brune 
413fde0ff24SPeter Brune    Logically Collective on SNES
414c79ef259SPeter Brune 
415c79ef259SPeter Brune    Input Parameters:
416c79ef259SPeter Brune +  snes - the multigrid context
417c79ef259SPeter Brune -  n    - the number of smoothing steps
418c79ef259SPeter Brune 
419c79ef259SPeter Brune    Options Database Key:
420d28543b3SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
421c79ef259SPeter Brune 
422c79ef259SPeter Brune    Level: advanced
423c79ef259SPeter Brune 
424c79ef259SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
425c79ef259SPeter Brune 
426c79ef259SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
427c79ef259SPeter Brune @*/
428c79ef259SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
429c79ef259SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
430c79ef259SPeter Brune   PetscErrorCode ierr = 0;
431c79ef259SPeter Brune   PetscFunctionBegin;
432c79ef259SPeter Brune   fas->max_down_it = n;
433c79ef259SPeter Brune   if (fas->next) {
434c79ef259SPeter Brune     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
435c79ef259SPeter Brune   }
436c79ef259SPeter Brune   PetscFunctionReturn(0);
437c79ef259SPeter Brune }
438c79ef259SPeter Brune 
439c79ef259SPeter Brune #undef __FUNCT__
440421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
441c79ef259SPeter Brune /*@
442c79ef259SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
443c79ef259SPeter Brune    interpolation from l-1 to the lth level
444c79ef259SPeter Brune 
445c79ef259SPeter Brune    Input Parameters:
446c79ef259SPeter Brune +  snes      - the multigrid context
447c79ef259SPeter Brune .  mat       - the interpolation operator
448c79ef259SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
449c79ef259SPeter Brune 
450c79ef259SPeter Brune    Level: advanced
451c79ef259SPeter Brune 
452c79ef259SPeter Brune    Notes:
453c79ef259SPeter Brune           Usually this is the same matrix used also to set the restriction
454c79ef259SPeter Brune     for the same level.
455c79ef259SPeter Brune 
456c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
457c79ef259SPeter Brune     out from the matrix size which one it is.
458c79ef259SPeter Brune 
459c79ef259SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
460c79ef259SPeter Brune 
461bccf9bb3SJed Brown .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
462c79ef259SPeter Brune @*/
463421d9b32SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
464421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
465d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
466bccf9bb3SJed Brown   PetscErrorCode ierr;
467421d9b32SPeter Brune 
468421d9b32SPeter Brune   PetscFunctionBegin;
469*2e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D interpolation from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
470421d9b32SPeter Brune   /* get to the correct level */
471d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
472421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
473421d9b32SPeter Brune   }
474*2e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
475bccf9bb3SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
476421d9b32SPeter Brune   fas->interpolate = mat;
477421d9b32SPeter Brune   PetscFunctionReturn(0);
478421d9b32SPeter Brune }
479421d9b32SPeter Brune 
480421d9b32SPeter Brune #undef __FUNCT__
481421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
482c79ef259SPeter Brune /*@
483c79ef259SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
484c79ef259SPeter Brune    from level l to l-1.
485c79ef259SPeter Brune 
486c79ef259SPeter Brune    Input Parameters:
487c79ef259SPeter Brune +  snes  - the multigrid context
488c79ef259SPeter Brune .  mat   - the restriction matrix
489c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
490c79ef259SPeter Brune 
491c79ef259SPeter Brune    Level: advanced
492c79ef259SPeter Brune 
493c79ef259SPeter Brune    Notes:
494c79ef259SPeter Brune           Usually this is the same matrix used also to set the interpolation
495c79ef259SPeter Brune     for the same level.
496c79ef259SPeter Brune 
497c79ef259SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
498c79ef259SPeter Brune     out from the matrix size which one it is.
499c79ef259SPeter Brune 
500fde0ff24SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
501c79ef259SPeter Brune     is used.
502c79ef259SPeter Brune 
503c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
504c79ef259SPeter Brune 
505c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
506c79ef259SPeter Brune @*/
507421d9b32SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
508421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
509d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
510bccf9bb3SJed Brown   PetscErrorCode ierr;
511421d9b32SPeter Brune 
512421d9b32SPeter Brune   PetscFunctionBegin;
513*2e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D restriction from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
514421d9b32SPeter Brune   /* get to the correct level */
515d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
516421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
517421d9b32SPeter Brune   }
518*2e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
519bccf9bb3SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
520421d9b32SPeter Brune   fas->restrct = mat;
521421d9b32SPeter Brune   PetscFunctionReturn(0);
522421d9b32SPeter Brune }
523421d9b32SPeter Brune 
524421d9b32SPeter Brune #undef __FUNCT__
525421d9b32SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
526c79ef259SPeter Brune /*@
527c79ef259SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
528c79ef259SPeter Brune    operator from level l to l-1.
529c79ef259SPeter Brune 
530c79ef259SPeter Brune    Input Parameters:
531c79ef259SPeter Brune +  snes   - the multigrid context
532c79ef259SPeter Brune .  rscale - the restriction scaling
533c79ef259SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
534c79ef259SPeter Brune 
535c79ef259SPeter Brune    Level: advanced
536c79ef259SPeter Brune 
537c79ef259SPeter Brune    Notes:
538c79ef259SPeter Brune          This is only used in the case that the injection is not set.
539c79ef259SPeter Brune 
540c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
541c79ef259SPeter Brune 
542c79ef259SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
543c79ef259SPeter Brune @*/
544421d9b32SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
545421d9b32SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
546d9ce41dcSBarry Smith   PetscInt top_level = fas->level,i;
547bccf9bb3SJed Brown   PetscErrorCode ierr;
548421d9b32SPeter Brune 
549421d9b32SPeter Brune   PetscFunctionBegin;
550*2e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D rscale from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
551421d9b32SPeter Brune   /* get to the correct level */
552d9ce41dcSBarry Smith   for (i = fas->level; i > level; i--) {
553421d9b32SPeter Brune     fas = (SNES_FAS *)fas->next->data;
554421d9b32SPeter Brune   }
555*2e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
556bccf9bb3SJed Brown   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
557421d9b32SPeter Brune   fas->rscale = rscale;
558421d9b32SPeter Brune   PetscFunctionReturn(0);
559421d9b32SPeter Brune }
560421d9b32SPeter Brune 
561efe1f98aSPeter Brune 
562efe1f98aSPeter Brune #undef __FUNCT__
563efe1f98aSPeter Brune #define __FUNCT__ "SNESFASSetInjection"
564c79ef259SPeter Brune /*@
565c79ef259SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
566c79ef259SPeter Brune    from level l to l-1.
567c79ef259SPeter Brune 
568c79ef259SPeter Brune    Input Parameters:
569c79ef259SPeter Brune +  snes  - the multigrid context
570c79ef259SPeter Brune .  mat   - the restriction matrix
571c79ef259SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
572c79ef259SPeter Brune 
573c79ef259SPeter Brune    Level: advanced
574c79ef259SPeter Brune 
575c79ef259SPeter Brune    Notes:
576c79ef259SPeter Brune          If you do not set this, the restriction and rscale is used to
577c79ef259SPeter Brune    project the solution instead.
578c79ef259SPeter Brune 
579c79ef259SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
580c79ef259SPeter Brune 
581c79ef259SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
582c79ef259SPeter Brune @*/
583efe1f98aSPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
584efe1f98aSPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
585efe1f98aSPeter Brune   PetscInt top_level = fas->level,i;
586bccf9bb3SJed Brown   PetscErrorCode ierr;
587efe1f98aSPeter Brune 
588efe1f98aSPeter Brune   PetscFunctionBegin;
589*2e8ce248SJed Brown   if (level > top_level) SETERRQ3(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Cannot set level %D injection from level %D of SNESFAS with %D levels total",level,top_level,fas->levels);
590efe1f98aSPeter Brune   /* get to the correct level */
591efe1f98aSPeter Brune   for (i = fas->level; i > level; i--) {
592efe1f98aSPeter Brune     fas = (SNES_FAS *)fas->next->data;
593efe1f98aSPeter Brune   }
594*2e8ce248SJed Brown   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
595bccf9bb3SJed Brown   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
596efe1f98aSPeter Brune   fas->inject = mat;
597efe1f98aSPeter Brune   PetscFunctionReturn(0);
598efe1f98aSPeter Brune }
599efe1f98aSPeter Brune 
600421d9b32SPeter Brune #undef __FUNCT__
601421d9b32SPeter Brune #define __FUNCT__ "SNESReset_FAS"
602421d9b32SPeter Brune PetscErrorCode SNESReset_FAS(SNES snes)
603421d9b32SPeter Brune {
60477df8cc4SPeter Brune   PetscErrorCode ierr = 0;
605421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
606421d9b32SPeter Brune 
607421d9b32SPeter Brune   PetscFunctionBegin;
608bccf9bb3SJed Brown   ierr = SNESDestroy(&fas->upsmooth);CHKERRQ(ierr);
609bccf9bb3SJed Brown   ierr = SNESDestroy(&fas->downsmooth);CHKERRQ(ierr);
6103dccd265SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
611bccf9bb3SJed Brown   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
612bccf9bb3SJed Brown   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
613bccf9bb3SJed Brown   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
6143dccd265SPeter Brune 
6153dccd265SPeter Brune   if (fas->upsmooth)   ierr = SNESReset(fas->upsmooth);CHKERRQ(ierr);
6163dccd265SPeter Brune   if (fas->downsmooth) ierr = SNESReset(fas->downsmooth);CHKERRQ(ierr);
617742fe5e2SPeter Brune   if (fas->next)       ierr = SNESReset(fas->next);CHKERRQ(ierr);
618421d9b32SPeter Brune   PetscFunctionReturn(0);
619421d9b32SPeter Brune }
620421d9b32SPeter Brune 
621421d9b32SPeter Brune #undef __FUNCT__
622421d9b32SPeter Brune #define __FUNCT__ "SNESDestroy_FAS"
623421d9b32SPeter Brune PetscErrorCode SNESDestroy_FAS(SNES snes)
624421d9b32SPeter Brune {
625421d9b32SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
626742fe5e2SPeter Brune   PetscErrorCode ierr = 0;
627421d9b32SPeter Brune 
628421d9b32SPeter Brune   PetscFunctionBegin;
629421d9b32SPeter Brune   /* recursively resets and then destroys */
63079d9a41aSPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
631421d9b32SPeter Brune   if (fas->next)         ierr = SNESDestroy(&fas->next);CHKERRQ(ierr);
632421d9b32SPeter Brune   ierr = PetscFree(fas);CHKERRQ(ierr);
633ee1fd11aSPeter Brune 
634421d9b32SPeter Brune   PetscFunctionReturn(0);
635421d9b32SPeter Brune }
636421d9b32SPeter Brune 
637421d9b32SPeter Brune #undef __FUNCT__
638421d9b32SPeter Brune #define __FUNCT__ "SNESSetUp_FAS"
639421d9b32SPeter Brune PetscErrorCode SNESSetUp_FAS(SNES snes)
640421d9b32SPeter Brune {
64148bfdf8aSPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
642421d9b32SPeter Brune   PetscErrorCode ierr;
643efe1f98aSPeter Brune   VecScatter     injscatter;
644d1adcc6fSPeter Brune   PetscInt       dm_levels;
6453dccd265SPeter Brune   Vec            vec_sol, vec_func, vec_sol_update, vec_rhs; /* preserve these if they're set through the reset */
646d1adcc6fSPeter Brune 
647421d9b32SPeter Brune   PetscFunctionBegin;
648eff52c0eSPeter Brune 
649cc05f883SPeter Brune   if (fas->usedmfornumberoflevels && (fas->level == fas->levels - 1)) {
650d1adcc6fSPeter Brune     ierr = DMGetRefineLevel(snes->dm,&dm_levels);CHKERRQ(ierr);
651d1adcc6fSPeter Brune     dm_levels++;
652cc05f883SPeter Brune     if (dm_levels > fas->levels) {
6533dccd265SPeter Brune 
654*2e8ce248SJed Brown       /* we don't want the solution and func vectors to be destroyed in the SNESReset when it's called in SNESFASSetLevels_FAS*/
6553dccd265SPeter Brune       vec_sol = snes->vec_sol;
6563dccd265SPeter Brune       vec_func = snes->vec_func;
6573dccd265SPeter Brune       vec_sol_update = snes->vec_sol_update;
6583dccd265SPeter Brune       vec_rhs = snes->vec_rhs;
6593dccd265SPeter Brune       snes->vec_sol = PETSC_NULL;
6603dccd265SPeter Brune       snes->vec_func = PETSC_NULL;
6613dccd265SPeter Brune       snes->vec_sol_update = PETSC_NULL;
6623dccd265SPeter Brune       snes->vec_rhs = PETSC_NULL;
6633dccd265SPeter Brune 
6643dccd265SPeter Brune       /* reset the number of levels */
665d1adcc6fSPeter Brune       ierr = SNESFASSetLevels(snes,dm_levels,PETSC_NULL);CHKERRQ(ierr);
666cc05f883SPeter Brune       ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
6673dccd265SPeter Brune 
6683dccd265SPeter Brune       snes->vec_sol = vec_sol;
6693dccd265SPeter Brune       snes->vec_func = vec_func;
6703dccd265SPeter Brune       snes->vec_rhs = vec_rhs;
6713dccd265SPeter Brune       snes->vec_sol_update = vec_sol_update;
672d1adcc6fSPeter Brune     }
673d1adcc6fSPeter Brune   }
674d1adcc6fSPeter Brune 
6753dccd265SPeter Brune   if (fas->level != fas->levels - 1) snes->gridsequence = 0; /* no grid sequencing inside the multigrid hierarchy! */
6763dccd265SPeter Brune 
67707144faaSPeter Brune   if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
678fde0ff24SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
67907144faaSPeter Brune   } else {
680ddebd997SPeter Brune     ierr = SNESDefaultGetWork(snes, 4);CHKERRQ(ierr); /* work vectors used for intergrid transfers */
68107144faaSPeter Brune   }
682cc05f883SPeter Brune 
68379d9a41aSPeter Brune   if (snes->dm) {
68479d9a41aSPeter Brune     /* construct EVERYTHING from the DM -- including the progressive set of smoothers */
68579d9a41aSPeter Brune     if (fas->next) {
68679d9a41aSPeter Brune       /* for now -- assume the DM and the evaluation functions have been set externally */
68779d9a41aSPeter Brune       if (!fas->next->dm) {
68879d9a41aSPeter Brune         ierr = DMCoarsen(snes->dm, ((PetscObject)fas->next)->comm, &fas->next->dm);CHKERRQ(ierr);
68979d9a41aSPeter Brune         ierr = SNESSetDM(fas->next, fas->next->dm);CHKERRQ(ierr);
69079d9a41aSPeter Brune       }
69179d9a41aSPeter Brune       /* set the interpolation and restriction from the DM */
69279d9a41aSPeter Brune       if (!fas->interpolate) {
69379d9a41aSPeter Brune         ierr = DMCreateInterpolation(fas->next->dm, snes->dm, &fas->interpolate, &fas->rscale);CHKERRQ(ierr);
694bccf9bb3SJed Brown         if (!fas->restrct) {
695bccf9bb3SJed Brown           ierr = PetscObjectReference((PetscObject)fas->interpolate);CHKERRQ(ierr);
69679d9a41aSPeter Brune           fas->restrct = fas->interpolate;
69779d9a41aSPeter Brune         }
698bccf9bb3SJed Brown       }
69979d9a41aSPeter Brune       /* set the injection from the DM */
70079d9a41aSPeter Brune       if (!fas->inject) {
70179d9a41aSPeter Brune         ierr = DMCreateInjection(fas->next->dm, snes->dm, &injscatter);CHKERRQ(ierr);
70279d9a41aSPeter Brune         ierr = MatCreateScatter(((PetscObject)snes)->comm, injscatter, &fas->inject);CHKERRQ(ierr);
70379d9a41aSPeter Brune         ierr = VecScatterDestroy(&injscatter);CHKERRQ(ierr);
70479d9a41aSPeter Brune       }
70579d9a41aSPeter Brune     }
70679d9a41aSPeter Brune     /* set the DMs of the pre and post-smoothers here */
70779d9a41aSPeter Brune     if (fas->upsmooth)  {ierr = SNESSetDM(fas->upsmooth,   snes->dm);CHKERRQ(ierr);}
70879d9a41aSPeter Brune     if (fas->downsmooth){ierr = SNESSetDM(fas->downsmooth, snes->dm);CHKERRQ(ierr);}
70979d9a41aSPeter Brune   }
71079d9a41aSPeter Brune   /*pass the smoother, function, and jacobian up to the next level if it's not user set already */
71179d9a41aSPeter Brune 
71279d9a41aSPeter Brune  if (fas->next) {
71379d9a41aSPeter Brune     if (fas->galerkin) {
71479d9a41aSPeter Brune       ierr = SNESSetFunction(fas->next, PETSC_NULL, SNESFASGalerkinDefaultFunction, fas->next);CHKERRQ(ierr);
71579d9a41aSPeter Brune     } else {
71679d9a41aSPeter Brune       if (snes->ops->computefunction && !fas->next->ops->computefunction) {
71779d9a41aSPeter Brune         ierr = SNESSetFunction(fas->next, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
71879d9a41aSPeter Brune       }
71979d9a41aSPeter Brune       if (snes->ops->computejacobian && !fas->next->ops->computejacobian) {
72079d9a41aSPeter Brune         ierr = SNESSetJacobian(fas->next, fas->next->jacobian, fas->next->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
72179d9a41aSPeter Brune       }
72279d9a41aSPeter Brune       if (snes->ops->computegs && !fas->next->ops->computegs) {
72379d9a41aSPeter Brune         ierr = SNESSetGS(fas->next, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
72479d9a41aSPeter Brune       }
72579d9a41aSPeter Brune     }
72679d9a41aSPeter Brune   }
72779d9a41aSPeter Brune 
72879d9a41aSPeter Brune   if (fas->next) {
72979d9a41aSPeter Brune     /* gotta set up the solution vector for this to work */
73079d9a41aSPeter Brune     if (!fas->next->vec_sol) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_sol);CHKERRQ(ierr);}
73179d9a41aSPeter Brune     if (!fas->next->vec_rhs) {ierr = VecDuplicate(fas->rscale, &fas->next->vec_rhs);CHKERRQ(ierr);}
73279d9a41aSPeter Brune     ierr = SNESSetUp(fas->next);CHKERRQ(ierr);
73379d9a41aSPeter Brune   }
73479d9a41aSPeter Brune 
73548bfdf8aSPeter Brune   /* setup the pre and post smoothers and set their function, jacobian, and gs evaluation routines if the user has neglected this */
73648bfdf8aSPeter Brune   if (fas->upsmooth) {
73748bfdf8aSPeter Brune     if (snes->ops->computefunction && !fas->upsmooth->ops->computefunction) {
73848bfdf8aSPeter Brune       ierr = SNESSetFunction(fas->upsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
73948bfdf8aSPeter Brune     }
74048bfdf8aSPeter Brune     if (snes->ops->computejacobian && !fas->upsmooth->ops->computejacobian) {
74148bfdf8aSPeter Brune       ierr = SNESSetJacobian(fas->upsmooth, fas->upsmooth->jacobian, fas->upsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
74248bfdf8aSPeter Brune     }
74348bfdf8aSPeter Brune    if (snes->ops->computegs && !fas->upsmooth->ops->computegs) {
74448bfdf8aSPeter Brune       ierr = SNESSetGS(fas->upsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
74548bfdf8aSPeter Brune     }
746fde0ff24SPeter Brune    ierr = SNESSetFromOptions(fas->upsmooth);CHKERRQ(ierr);
74748bfdf8aSPeter Brune   }
74848bfdf8aSPeter Brune   if (fas->downsmooth) {
74948bfdf8aSPeter Brune     if (snes->ops->computefunction && !fas->downsmooth->ops->computefunction) {
75048bfdf8aSPeter Brune       ierr = SNESSetFunction(fas->downsmooth, PETSC_NULL, snes->ops->computefunction, snes->funP);CHKERRQ(ierr);
75148bfdf8aSPeter Brune     }
75248bfdf8aSPeter Brune     if (snes->ops->computejacobian && !fas->downsmooth->ops->computejacobian) {
75348bfdf8aSPeter Brune       ierr = SNESSetJacobian(fas->downsmooth, fas->downsmooth->jacobian, fas->downsmooth->jacobian_pre, snes->ops->computejacobian, snes->jacP);CHKERRQ(ierr);
75448bfdf8aSPeter Brune     }
75548bfdf8aSPeter Brune    if (snes->ops->computegs && !fas->downsmooth->ops->computegs) {
75648bfdf8aSPeter Brune      ierr = SNESSetGS(fas->downsmooth, snes->ops->computegs, snes->gsP);CHKERRQ(ierr);
75748bfdf8aSPeter Brune     }
758fde0ff24SPeter Brune     ierr = SNESSetFromOptions(fas->downsmooth);CHKERRQ(ierr);
7591a266240SBarry Smith   }
760cc05f883SPeter Brune 
7616273346dSPeter Brune   /* setup FAS work vectors */
7626273346dSPeter Brune   if (fas->galerkin) {
7636273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Xg);CHKERRQ(ierr);
7646273346dSPeter Brune     ierr = VecDuplicate(snes->vec_sol, &fas->Fg);CHKERRQ(ierr);
7656273346dSPeter Brune   }
7666273346dSPeter Brune 
76779d9a41aSPeter Brune 
768fa9694d7SPeter Brune   /* got to set them all up at once */
769421d9b32SPeter Brune   PetscFunctionReturn(0);
770421d9b32SPeter Brune }
771421d9b32SPeter Brune 
772421d9b32SPeter Brune #undef __FUNCT__
773421d9b32SPeter Brune #define __FUNCT__ "SNESSetFromOptions_FAS"
774421d9b32SPeter Brune PetscErrorCode SNESSetFromOptions_FAS(SNES snes)
775421d9b32SPeter Brune {
776ee78dd50SPeter Brune   SNES_FAS       *fas = (SNES_FAS *) snes->data;
777ee78dd50SPeter Brune   PetscInt       levels = 1;
778fde0ff24SPeter Brune   PetscBool      flg, smoothflg, smoothupflg, smoothdownflg, smoothcoarseflg = PETSC_FALSE, monflg;
779421d9b32SPeter Brune   PetscErrorCode ierr;
780ee78dd50SPeter Brune   char           monfilename[PETSC_MAX_PATH_LEN];
78107144faaSPeter Brune   SNESFASType    fastype;
782fde0ff24SPeter Brune   const char     *optionsprefix;
783fde0ff24SPeter Brune   const char     *prefix;
784421d9b32SPeter Brune 
785421d9b32SPeter Brune   PetscFunctionBegin;
786c90fad12SPeter Brune   ierr = PetscOptionsHead("SNESFAS Options-----------------------------------");CHKERRQ(ierr);
787ee78dd50SPeter Brune 
788ee78dd50SPeter Brune   /* number of levels -- only process on the finest level */
789ee78dd50SPeter Brune   if (fas->levels - 1 == fas->level) {
790ee78dd50SPeter Brune     ierr = PetscOptionsInt("-snes_fas_levels", "Number of Levels", "SNESFASSetLevels", levels, &levels, &flg);CHKERRQ(ierr);
791c732cbdbSBarry Smith     if (!flg && snes->dm) {
792c732cbdbSBarry Smith       ierr = DMGetRefineLevel(snes->dm,&levels);CHKERRQ(ierr);
793c732cbdbSBarry Smith       levels++;
794d1adcc6fSPeter Brune       fas->usedmfornumberoflevels = PETSC_TRUE;
795c732cbdbSBarry Smith     }
796ee78dd50SPeter Brune     ierr = SNESFASSetLevels(snes, levels, PETSC_NULL);CHKERRQ(ierr);
797ee78dd50SPeter Brune   }
79807144faaSPeter Brune   fastype = fas->fastype;
79907144faaSPeter Brune   ierr = PetscOptionsEnum("-snes_fas_type","FAS correction type","SNESFASSetType",SNESFASTypes,(PetscEnum)fastype,(PetscEnum*)&fastype,&flg);CHKERRQ(ierr);
80007144faaSPeter Brune   if (flg) {
80107144faaSPeter Brune     ierr = SNESFASSetType(snes, fastype);CHKERRQ(ierr);
80207144faaSPeter Brune   }
803ee78dd50SPeter Brune 
804fde0ff24SPeter Brune   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
805fde0ff24SPeter Brune 
806fde0ff24SPeter Brune   /* smoother setup options */
807fde0ff24SPeter Brune   ierr = PetscOptionsHasName(optionsprefix, "-fas_snes_type", &smoothflg);CHKERRQ(ierr);
808fde0ff24SPeter Brune   ierr = PetscOptionsHasName(optionsprefix, "-fas_up_snes_type", &smoothupflg);CHKERRQ(ierr);
809fde0ff24SPeter Brune   ierr = PetscOptionsHasName(optionsprefix, "-fas_down_snes_type", &smoothdownflg);CHKERRQ(ierr);
810fde0ff24SPeter Brune   if (fas->level == 0) {
811fde0ff24SPeter Brune     ierr = PetscOptionsHasName(optionsprefix, "-fas_coarse_snes_type", &smoothcoarseflg);CHKERRQ(ierr);
812fde0ff24SPeter Brune   }
813ee78dd50SPeter Brune   /* options for the number of preconditioning cycles and cycle type */
814fde0ff24SPeter Brune 
815d1adcc6fSPeter Brune   ierr = PetscOptionsInt("-snes_fas_cycles","Number of cycles","SNESFASSetCycles",fas->n_cycles,&fas->n_cycles,&flg);CHKERRQ(ierr);
816ee78dd50SPeter Brune 
8176273346dSPeter Brune   ierr = PetscOptionsBool("-snes_fas_galerkin", "Form coarse problems with Galerkin","SNESFAS",fas->galerkin,&fas->galerkin,&flg);CHKERRQ(ierr);
8186273346dSPeter Brune 
819646217ecSPeter Brune   ierr = PetscOptionsString("-snes_fas_monitor","Monitor FAS progress","SNESMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&monflg);CHKERRQ(ierr);
820ee78dd50SPeter Brune 
821ee78dd50SPeter Brune 
822421d9b32SPeter Brune   ierr = PetscOptionsTail();CHKERRQ(ierr);
8238cc86e31SPeter Brune   /* setup from the determined types if there is no pointwise procedure or smoother defined */
824eff52c0eSPeter Brune 
825fde0ff24SPeter Brune   if ((!fas->downsmooth) && smoothcoarseflg) {
826fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
827fde0ff24SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
828fde0ff24SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
829fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_coarse_");CHKERRQ(ierr);
830fde0ff24SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
831fde0ff24SPeter Brune   }
832fde0ff24SPeter Brune 
833fde0ff24SPeter Brune   if ((!fas->downsmooth) && smoothdownflg) {
8348cc86e31SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
8358cc86e31SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
8368cc86e31SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
837fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_down_");CHKERRQ(ierr);
8388cc86e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
8398cc86e31SPeter Brune   }
8408cc86e31SPeter Brune 
841fde0ff24SPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && smoothupflg) {
84267339d5cSBarry Smith     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
843ee78dd50SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
84467339d5cSBarry Smith     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
845fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_up_");CHKERRQ(ierr);
846293a7e31SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
847ee78dd50SPeter Brune   }
848fde0ff24SPeter Brune 
849fde0ff24SPeter Brune   if ((!fas->downsmooth) && smoothflg) {
850fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
851fde0ff24SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->downsmooth);CHKERRQ(ierr);
852fde0ff24SPeter Brune     ierr = SNESSetOptionsPrefix(fas->downsmooth,prefix);CHKERRQ(ierr);
853fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->downsmooth,"fas_");CHKERRQ(ierr);
854fde0ff24SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->downsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
855fde0ff24SPeter Brune   }
856fde0ff24SPeter Brune 
857fde0ff24SPeter Brune   if ((!fas->upsmooth) && (fas->level != 0) && smoothflg) {
858fde0ff24SPeter Brune     ierr = SNESGetOptionsPrefix(snes,&prefix);CHKERRQ(ierr);
859fde0ff24SPeter Brune     ierr = SNESCreate(((PetscObject)snes)->comm, &fas->upsmooth);CHKERRQ(ierr);
860fde0ff24SPeter Brune     ierr = SNESSetOptionsPrefix(fas->upsmooth,prefix);CHKERRQ(ierr);
861fde0ff24SPeter Brune     ierr = SNESAppendOptionsPrefix(fas->upsmooth,"fas_");CHKERRQ(ierr);
862fde0ff24SPeter Brune     ierr = PetscObjectIncrementTabLevel((PetscObject)fas->upsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
863fde0ff24SPeter Brune   }
864fde0ff24SPeter Brune 
8651a266240SBarry Smith   if (fas->upsmooth) {
866fde0ff24SPeter Brune     ierr = SNESSetTolerances(fas->upsmooth, fas->upsmooth->abstol, fas->upsmooth->rtol, fas->upsmooth->xtol, 1, 1000);CHKERRQ(ierr);
8671a266240SBarry Smith   }
8681a266240SBarry Smith 
8691a266240SBarry Smith   if (fas->downsmooth) {
870fde0ff24SPeter Brune     ierr = SNESSetTolerances(fas->downsmooth, fas->downsmooth->abstol, fas->downsmooth->rtol, fas->downsmooth->xtol, 1, 1000);CHKERRQ(ierr);
8711a266240SBarry Smith   }
872ee78dd50SPeter Brune 
873742fe5e2SPeter Brune   if (fas->level != fas->levels - 1) {
874fde0ff24SPeter Brune     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->xtol, fas->n_cycles, 1000);CHKERRQ(ierr);
875742fe5e2SPeter Brune   }
876742fe5e2SPeter Brune 
877ce11c761SPeter Brune   /* control the simple Richardson smoother that is default if there's no upsmooth or downsmooth */
878ce11c761SPeter Brune   if (!fas->downsmooth) {
87979d9a41aSPeter Brune     ierr = PetscOptionsInt("-fas_down_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
88093dfaa4eSPeter Brune     ierr = PetscOptionsInt("-fas_snes_max_it","Down smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
881ce11c761SPeter Brune     if (fas->level == 0) {
88279d9a41aSPeter Brune       ierr = PetscOptionsInt("-fas_coarse_snes_max_it","Coarse smooth iterations","SNESFASSetCycles",fas->max_down_it,&fas->max_down_it,&flg);CHKERRQ(ierr);
883ce11c761SPeter Brune     }
884ce11c761SPeter Brune   }
885ce11c761SPeter Brune 
886ce11c761SPeter Brune   if (!fas->upsmooth) {
88779d9a41aSPeter Brune     ierr = PetscOptionsInt("-fas_up_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
88893dfaa4eSPeter Brune     ierr = PetscOptionsInt("-fas_snes_max_it","Upsmooth iterations","SNESFASSetCycles",fas->max_up_it,&fas->max_up_it,&flg);CHKERRQ(ierr);
889ce11c761SPeter Brune   }
890ce11c761SPeter Brune 
891ee78dd50SPeter Brune   if (monflg) {
892646217ecSPeter Brune     fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)snes)->comm);CHKERRQ(ierr);
893794bee33SPeter Brune     /* set the monitors for the upsmoother and downsmoother */
8942f7ea302SPeter Brune     ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
895742fe5e2SPeter Brune     ierr = SNESMonitorSet(snes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
896293a7e31SPeter Brune     if (fas->upsmooth)   ierr = SNESMonitorSet(fas->upsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
897293a7e31SPeter Brune     if (fas->downsmooth) ierr = SNESMonitorSet(fas->downsmooth,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
898d28543b3SPeter Brune   } else {
899d28543b3SPeter Brune     /* unset the monitors on the coarse levels */
900d28543b3SPeter Brune     if (fas->level != fas->levels - 1) {
901d28543b3SPeter Brune       ierr = SNESMonitorCancel(snes);CHKERRQ(ierr);
902d28543b3SPeter Brune     }
903ee78dd50SPeter Brune   }
904ee78dd50SPeter Brune 
905ee78dd50SPeter Brune   /* recursive option setting for the smoothers */
906d28543b3SPeter Brune   if (fas->next) {ierr = SNESSetFromOptions(fas->next);CHKERRQ(ierr);}
907421d9b32SPeter Brune   PetscFunctionReturn(0);
908421d9b32SPeter Brune }
909421d9b32SPeter Brune 
910421d9b32SPeter Brune #undef __FUNCT__
911421d9b32SPeter Brune #define __FUNCT__ "SNESView_FAS"
912421d9b32SPeter Brune PetscErrorCode SNESView_FAS(SNES snes, PetscViewer viewer)
913421d9b32SPeter Brune {
914421d9b32SPeter Brune   SNES_FAS   *fas = (SNES_FAS *) snes->data;
915421d9b32SPeter Brune   PetscBool      iascii;
916421d9b32SPeter Brune   PetscErrorCode ierr;
917421d9b32SPeter Brune 
918421d9b32SPeter Brune   PetscFunctionBegin;
919421d9b32SPeter Brune   ierr = PetscTypeCompare((PetscObject) viewer, PETSCVIEWERASCII, &iascii);CHKERRQ(ierr);
920421d9b32SPeter Brune   if (iascii) {
921*2e8ce248SJed Brown     ierr = PetscViewerASCIIPrintf(viewer, "FAS, levels = %D\n",  fas->levels);CHKERRQ(ierr);
922421d9b32SPeter Brune     ierr = PetscViewerASCIIPushTab(viewer);
923421d9b32SPeter Brune     ierr = PetscViewerASCIIPrintf(viewer, "level: %d\n",  fas->level);CHKERRQ(ierr);
924ee78dd50SPeter Brune     if (fas->upsmooth) {
92539a4b097SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "up-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
926421d9b32SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);
927ee78dd50SPeter Brune       ierr = SNESView(fas->upsmooth, viewer);CHKERRQ(ierr);
928421d9b32SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);
929421d9b32SPeter Brune     } else {
930ee78dd50SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "no up-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
931421d9b32SPeter Brune     }
932ee78dd50SPeter Brune     if (fas->downsmooth) {
93339a4b097SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "down-smoother on level %D:\n",  fas->level);CHKERRQ(ierr);
934421d9b32SPeter Brune       ierr = PetscViewerASCIIPushTab(viewer);
935ee78dd50SPeter Brune       ierr = SNESView(fas->downsmooth, viewer);CHKERRQ(ierr);
936421d9b32SPeter Brune       ierr = PetscViewerASCIIPopTab(viewer);
937421d9b32SPeter Brune     } else {
938ee78dd50SPeter Brune       ierr = PetscViewerASCIIPrintf(viewer, "no down-smoother on level %D\n",  fas->level);CHKERRQ(ierr);
939421d9b32SPeter Brune     }
940421d9b32SPeter Brune     ierr = PetscViewerASCIIPopTab(viewer);
941421d9b32SPeter Brune   } else {
942421d9b32SPeter Brune     SETERRQ1(((PetscObject)snes)->comm,PETSC_ERR_SUP,"Viewer type %s not supported for SNESFAS",((PetscObject)viewer)->type_name);
943421d9b32SPeter Brune   }
944421d9b32SPeter Brune   PetscFunctionReturn(0);
945421d9b32SPeter Brune }
946421d9b32SPeter Brune 
947421d9b32SPeter Brune #undef __FUNCT__
94839bd7f45SPeter Brune #define __FUNCT__ "FASDownSmooth"
94939bd7f45SPeter Brune /*
95039bd7f45SPeter Brune Defines the action of the downsmoother
95139bd7f45SPeter Brune  */
95239bd7f45SPeter Brune PetscErrorCode FASDownSmooth(SNES snes, Vec B, Vec X, Vec F){
95339bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
954fde0ff24SPeter Brune   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
955742fe5e2SPeter Brune   SNESConvergedReason reason;
95639bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
9574b32a720SPeter Brune   Vec                 G, W, Y, FPC;
958fde0ff24SPeter Brune   PetscBool           lssuccess;
959fde0ff24SPeter Brune   PetscInt            k;
960421d9b32SPeter Brune   PetscFunctionBegin;
961fde0ff24SPeter Brune   G = snes->work[1];
962fde0ff24SPeter Brune   W = snes->work[2];
963fde0ff24SPeter Brune   Y = snes->work[3];
964d1adcc6fSPeter Brune   if (fas->downsmooth) {
965d1adcc6fSPeter Brune     ierr = SNESSolve(fas->downsmooth, B, X);CHKERRQ(ierr);
966742fe5e2SPeter Brune     /* check convergence reason for the smoother */
967742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->downsmooth,&reason);CHKERRQ(ierr);
968742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
969742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
970742fe5e2SPeter Brune       PetscFunctionReturn(0);
971742fe5e2SPeter Brune     }
9724b32a720SPeter Brune     ierr = SNESGetFunction(fas->downsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
9734b32a720SPeter Brune     ierr = VecCopy(FPC, F);CHKERRQ(ierr);
974fde0ff24SPeter Brune   } else {
975fde0ff24SPeter Brune     /* basic richardson smoother */
976fde0ff24SPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
977794bee33SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
978794bee33SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
979fde0ff24SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
980fde0ff24SPeter Brune       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
981fde0ff24SPeter Brune       if (!lssuccess) {
982fde0ff24SPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
983fde0ff24SPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
9842f7ea302SPeter Brune           PetscFunctionReturn(0);
9852f7ea302SPeter Brune         }
986fe6f9142SPeter Brune       }
987fde0ff24SPeter Brune       ierr = VecCopy(W, X);CHKERRQ(ierr);
988fde0ff24SPeter Brune       ierr = VecCopy(G, F);CHKERRQ(ierr);
989fde0ff24SPeter Brune       fnorm = gnorm;
990fde0ff24SPeter Brune     }
991fde0ff24SPeter Brune   }
99239bd7f45SPeter Brune   PetscFunctionReturn(0);
99339bd7f45SPeter Brune }
99439bd7f45SPeter Brune 
99539bd7f45SPeter Brune 
99639bd7f45SPeter Brune #undef __FUNCT__
99739bd7f45SPeter Brune #define __FUNCT__ "FASUpSmooth"
99839bd7f45SPeter Brune /*
99907144faaSPeter Brune Defines the action of the upsmoother
100039bd7f45SPeter Brune  */
100139bd7f45SPeter Brune PetscErrorCode FASUpSmooth (SNES snes, Vec B, Vec X, Vec F) {
100239bd7f45SPeter Brune   PetscErrorCode      ierr = 0;
1003fde0ff24SPeter Brune   PetscReal           fnorm, gnorm, ynorm, xnorm = 0.0;
100439bd7f45SPeter Brune   SNESConvergedReason reason;
100539bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
10064b32a720SPeter Brune   Vec                 G, W, Y, FPC;
1007fde0ff24SPeter Brune   PetscBool           lssuccess;
1008fde0ff24SPeter Brune   PetscInt            k;
100939bd7f45SPeter Brune   PetscFunctionBegin;
1010fde0ff24SPeter Brune   G = snes->work[1];
1011fde0ff24SPeter Brune   W = snes->work[2];
1012fde0ff24SPeter Brune   Y = snes->work[3];
101339bd7f45SPeter Brune   if (fas->upsmooth) {
1014fde0ff24SPeter Brune     ierr = SNESSolve(fas->upsmooth, B, X);CHKERRQ(ierr);
101539bd7f45SPeter Brune     /* check convergence reason for the smoother */
1016fde0ff24SPeter Brune     ierr = SNESGetConvergedReason(fas->upsmooth,&reason);CHKERRQ(ierr);
101739bd7f45SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
101839bd7f45SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
101939bd7f45SPeter Brune       PetscFunctionReturn(0);
102039bd7f45SPeter Brune     }
10214b32a720SPeter Brune     ierr = SNESGetFunction(fas->upsmooth, &FPC, PETSC_NULL, PETSC_NULL);CHKERRQ(ierr);
10224b32a720SPeter Brune     ierr = VecCopy(FPC, F);CHKERRQ(ierr);
1023fde0ff24SPeter Brune   } else {
1024fde0ff24SPeter Brune     /* basic richardson smoother */
1025fde0ff24SPeter Brune     for (k = 0; k < fas->max_up_it; k++) {
102639bd7f45SPeter Brune       ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
102739bd7f45SPeter Brune       ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1028fde0ff24SPeter Brune       ierr = VecCopy(F, Y);CHKERRQ(ierr);
1029fde0ff24SPeter Brune       ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Y,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssuccess);CHKERRQ(ierr);
1030fde0ff24SPeter Brune       if (!lssuccess) {
1031fde0ff24SPeter Brune         if (++snes->numFailures >= snes->maxFailures) {
1032fde0ff24SPeter Brune           snes->reason = SNES_DIVERGED_LINE_SEARCH;
103339bd7f45SPeter Brune           PetscFunctionReturn(0);
103439bd7f45SPeter Brune         }
103539bd7f45SPeter Brune       }
1036fde0ff24SPeter Brune       ierr = VecCopy(W, X);CHKERRQ(ierr);
1037fde0ff24SPeter Brune       ierr = VecCopy(G, F);CHKERRQ(ierr);
1038fde0ff24SPeter Brune       fnorm = gnorm;
1039fde0ff24SPeter Brune     }
1040fde0ff24SPeter Brune   }
104139bd7f45SPeter Brune   PetscFunctionReturn(0);
104239bd7f45SPeter Brune }
104339bd7f45SPeter Brune 
104439bd7f45SPeter Brune #undef __FUNCT__
104539bd7f45SPeter Brune #define __FUNCT__ "FASCoarseCorrection"
104639bd7f45SPeter Brune /*
104739bd7f45SPeter Brune 
104839bd7f45SPeter Brune Performs the FAS coarse correction as:
104939bd7f45SPeter Brune 
105039bd7f45SPeter Brune fine problem: F(x) = 0
105139bd7f45SPeter Brune coarse problem: F^c(x) = b^c
105239bd7f45SPeter Brune 
105339bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
105439bd7f45SPeter Brune 
105539bd7f45SPeter Brune with correction:
105639bd7f45SPeter Brune 
105739bd7f45SPeter Brune 
105839bd7f45SPeter Brune 
105939bd7f45SPeter Brune  */
106039a4b097SPeter Brune PetscErrorCode FASCoarseCorrection(SNES snes, Vec X, Vec F, Vec X_new) {
106139bd7f45SPeter Brune   PetscErrorCode      ierr;
106239bd7f45SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c;
106339bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
106439bd7f45SPeter Brune   SNESConvergedReason reason;
106539bd7f45SPeter Brune   PetscFunctionBegin;
1066fa9694d7SPeter Brune   if (fas->next) {
1067c90fad12SPeter Brune     X_c  = fas->next->vec_sol;
1068293a7e31SPeter Brune     Xo_c = fas->next->work[0];
1069c90fad12SPeter Brune     F_c  = fas->next->vec_func;
1070742fe5e2SPeter Brune     B_c  = fas->next->vec_rhs;
1071efe1f98aSPeter Brune 
1072efe1f98aSPeter Brune     /* inject the solution */
1073efe1f98aSPeter Brune     if (fas->inject) {
1074a3cb90a9SPeter Brune       ierr = MatRestrict(fas->inject, X, Xo_c);CHKERRQ(ierr);
1075efe1f98aSPeter Brune     } else {
1076a3cb90a9SPeter Brune       ierr = MatRestrict(fas->restrct, X, Xo_c);CHKERRQ(ierr);
1077a3cb90a9SPeter Brune       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
1078efe1f98aSPeter Brune     }
1079293a7e31SPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
1080293a7e31SPeter Brune 
1081293a7e31SPeter Brune     /* restrict the defect */
1082293a7e31SPeter Brune     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
1083293a7e31SPeter Brune 
1084c90fad12SPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
1085ee78dd50SPeter Brune     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
1086c90fad12SPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
1087742fe5e2SPeter Brune 
1088293a7e31SPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
1089c90fad12SPeter Brune 
1090ee78dd50SPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
1091ee78dd50SPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
1092c90fad12SPeter Brune 
1093c90fad12SPeter Brune     /* recurse to the next level */
1094f5a6d4f9SBarry Smith     fas->next->vec_rhs = B_c;
1095742fe5e2SPeter Brune     /* ierr = FASCycle_Private(fas->next, X_c);CHKERRQ(ierr); */
1096742fe5e2SPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
1097742fe5e2SPeter Brune     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
1098742fe5e2SPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
1099742fe5e2SPeter Brune       snes->reason = SNES_DIVERGED_INNER;
1100742fe5e2SPeter Brune       PetscFunctionReturn(0);
1101742fe5e2SPeter Brune     }
1102742fe5e2SPeter Brune     /* fas->next->vec_rhs = PETSC_NULL; */
1103ee78dd50SPeter Brune 
1104fa9694d7SPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
1105fa9694d7SPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
110639bd7f45SPeter Brune     ierr = MatInterpolateAdd(fas->interpolate, X_c, X, X_new);CHKERRQ(ierr);
1107293a7e31SPeter Brune   }
110839bd7f45SPeter Brune   PetscFunctionReturn(0);
110939bd7f45SPeter Brune }
111039bd7f45SPeter Brune 
111139bd7f45SPeter Brune #undef __FUNCT__
111239bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Additive"
111339bd7f45SPeter Brune /*
111439bd7f45SPeter Brune 
111539bd7f45SPeter Brune The additive cycle looks like:
111639bd7f45SPeter Brune 
111707144faaSPeter Brune xhat = x
111807144faaSPeter Brune xhat = dS(x, b)
111907144faaSPeter Brune x = coarsecorrection(xhat, b_d)
112007144faaSPeter Brune x = x + nu*(xhat - x);
112139bd7f45SPeter Brune (optional) x = uS(x, b)
112239bd7f45SPeter Brune 
112339bd7f45SPeter Brune With the coarse RHS (defect correction) as below.
112439bd7f45SPeter Brune 
112539bd7f45SPeter Brune  */
112639bd7f45SPeter Brune PetscErrorCode FASCycle_Additive(SNES snes, Vec X) {
112707144faaSPeter Brune   Vec                 F, B, Xhat;
1128ddebd997SPeter Brune   Vec                 X_c, Xo_c, F_c, B_c, G, W;
112939bd7f45SPeter Brune   PetscErrorCode      ierr;
113007144faaSPeter Brune   SNES_FAS *          fas = (SNES_FAS *)snes->data;
113107144faaSPeter Brune   SNESConvergedReason reason;
1132ddebd997SPeter Brune   PetscReal           xnorm = 0., fnorm = 0., gnorm = 0., ynorm = 0.;
1133ddebd997SPeter Brune   PetscBool           lssucceed;
113439bd7f45SPeter Brune   PetscFunctionBegin;
113539bd7f45SPeter Brune 
113639bd7f45SPeter Brune   F = snes->vec_func;
113739bd7f45SPeter Brune   B = snes->vec_rhs;
1138fde0ff24SPeter Brune   Xhat = snes->work[3];
1139fde0ff24SPeter Brune   G    = snes->work[1];
1140fde0ff24SPeter Brune   W    = snes->work[2];
114107144faaSPeter Brune   ierr = VecCopy(X, Xhat);CHKERRQ(ierr);
114207144faaSPeter Brune   /* recurse first */
114307144faaSPeter Brune   if (fas->next) {
114407144faaSPeter Brune     ierr = SNESComputeFunction(snes, Xhat, F);CHKERRQ(ierr);
114539bd7f45SPeter Brune 
114607144faaSPeter Brune     X_c  = fas->next->vec_sol;
114707144faaSPeter Brune     Xo_c = fas->next->work[0];
114807144faaSPeter Brune     F_c  = fas->next->vec_func;
114907144faaSPeter Brune     B_c  = fas->next->vec_rhs;
115039bd7f45SPeter Brune 
115107144faaSPeter Brune     /* inject the solution */
115207144faaSPeter Brune     if (fas->inject) {
115307144faaSPeter Brune       ierr = MatRestrict(fas->inject, Xhat, Xo_c);CHKERRQ(ierr);
115407144faaSPeter Brune     } else {
115507144faaSPeter Brune       ierr = MatRestrict(fas->restrct, Xhat, Xo_c);CHKERRQ(ierr);
115607144faaSPeter Brune       ierr = VecPointwiseMult(Xo_c, fas->rscale, Xo_c);CHKERRQ(ierr);
115707144faaSPeter Brune     }
115807144faaSPeter Brune     ierr = VecScale(F, -1.0);CHKERRQ(ierr);
115907144faaSPeter Brune 
116007144faaSPeter Brune     /* restrict the defect */
116107144faaSPeter Brune     ierr = MatRestrict(fas->restrct, F, B_c);CHKERRQ(ierr);
116207144faaSPeter Brune 
116307144faaSPeter Brune     /* solve the coarse problem corresponding to F^c(x^c) = b^c = Rb + F^c(Rx) - RF(x) */
116407144faaSPeter Brune     fas->next->vec_rhs = PETSC_NULL;                                           /*unset the RHS to evaluate function instead of residual*/
116507144faaSPeter Brune     ierr = SNESComputeFunction(fas->next, Xo_c, F_c);CHKERRQ(ierr);
116607144faaSPeter Brune 
116707144faaSPeter Brune     ierr = VecAXPY(B_c, 1.0, F_c);CHKERRQ(ierr);                               /* add F_c(X) to the RHS */
116807144faaSPeter Brune 
116907144faaSPeter Brune     /* set initial guess of the coarse problem to the projected fine solution */
117007144faaSPeter Brune     ierr = VecCopy(Xo_c, X_c);CHKERRQ(ierr);
117107144faaSPeter Brune 
117207144faaSPeter Brune     /* recurse */
117307144faaSPeter Brune     fas->next->vec_rhs = B_c;
117407144faaSPeter Brune     ierr = SNESSolve(fas->next, B_c, X_c);CHKERRQ(ierr);
117507144faaSPeter Brune 
117607144faaSPeter Brune     /* smooth on this level */
117707144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
117807144faaSPeter Brune 
117907144faaSPeter Brune     ierr = SNESGetConvergedReason(fas->next,&reason);CHKERRQ(ierr);
118007144faaSPeter Brune     if (reason < 0 && reason != SNES_DIVERGED_MAX_IT) {
118107144faaSPeter Brune       snes->reason = SNES_DIVERGED_INNER;
118207144faaSPeter Brune       PetscFunctionReturn(0);
118307144faaSPeter Brune     }
118407144faaSPeter Brune 
118507144faaSPeter Brune     /* correct as x <- x + I(x^c - Rx)*/
118607144faaSPeter Brune     ierr = VecAXPY(X_c, -1.0, Xo_c);CHKERRQ(ierr);
1187ddebd997SPeter Brune     ierr = MatInterpolate(fas->interpolate, X_c, Xhat);CHKERRQ(ierr);
118807144faaSPeter Brune 
1189ddebd997SPeter Brune     /* additive correction of the coarse direction*/
1190ddebd997SPeter Brune     ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1191ddebd997SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr);
1192eb1825c3SPeter Brune     ierr = VecScale(Xhat, -1.0);CHKERRQ(ierr);
1193ddebd997SPeter Brune     ierr = (*snes->ops->linesearch)(snes,snes->lsP,X,F,Xhat,fnorm,xnorm,G,W,&ynorm,&gnorm,&lssucceed);CHKERRQ(ierr);
1194ddebd997SPeter Brune     ierr = VecCopy(W, X);CHKERRQ(ierr);
1195ddebd997SPeter Brune     ierr = VecCopy(G, F);CHKERRQ(ierr);
1196ddebd997SPeter Brune     fnorm = gnorm;
119707144faaSPeter Brune   } else {
119807144faaSPeter Brune     ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
119907144faaSPeter Brune   }
120039bd7f45SPeter Brune   PetscFunctionReturn(0);
120139bd7f45SPeter Brune }
120239bd7f45SPeter Brune 
120339bd7f45SPeter Brune #undef __FUNCT__
120439bd7f45SPeter Brune #define __FUNCT__ "FASCycle_Multiplicative"
120539bd7f45SPeter Brune /*
120639bd7f45SPeter Brune 
120739bd7f45SPeter Brune Defines the FAS cycle as:
120839bd7f45SPeter Brune 
120939bd7f45SPeter Brune fine problem: F(x) = 0
121039bd7f45SPeter Brune coarse problem: F^c(x) = b^c
121139bd7f45SPeter Brune 
121239bd7f45SPeter Brune b^c = F^c(I^c_fx^f - I^c_fF(x))
121339bd7f45SPeter Brune 
121439bd7f45SPeter Brune correction:
121539bd7f45SPeter Brune 
121639bd7f45SPeter Brune x = x + I(x^c - Rx)
121739bd7f45SPeter Brune 
121839bd7f45SPeter Brune  */
121939bd7f45SPeter Brune PetscErrorCode FASCycle_Multiplicative(SNES snes, Vec X) {
122039bd7f45SPeter Brune 
122139bd7f45SPeter Brune   PetscErrorCode      ierr;
122239bd7f45SPeter Brune   Vec                 F,B;
122339bd7f45SPeter Brune   SNES_FAS            *fas = (SNES_FAS *)snes->data;
122439bd7f45SPeter Brune 
122539bd7f45SPeter Brune   PetscFunctionBegin;
122639bd7f45SPeter Brune   F = snes->vec_func;
122739bd7f45SPeter Brune   B = snes->vec_rhs;
122839bd7f45SPeter Brune   /* pre-smooth -- just update using the pre-smoother */
122939bd7f45SPeter Brune   ierr = FASDownSmooth(snes, B, X, F);CHKERRQ(ierr);
123039bd7f45SPeter Brune 
123139a4b097SPeter Brune   ierr = FASCoarseCorrection(snes, X, F, X);CHKERRQ(ierr);
123239bd7f45SPeter Brune 
1233c90fad12SPeter Brune   if (fas->level != 0) {
123439bd7f45SPeter Brune     ierr = FASUpSmooth(snes, B, X, F);CHKERRQ(ierr);
1235fe6f9142SPeter Brune   }
1236fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes, X, F);CHKERRQ(ierr);
1237fa9694d7SPeter Brune 
1238fa9694d7SPeter Brune   PetscFunctionReturn(0);
1239421d9b32SPeter Brune }
1240421d9b32SPeter Brune 
1241421d9b32SPeter Brune #undef __FUNCT__
1242421d9b32SPeter Brune #define __FUNCT__ "SNESSolve_FAS"
1243421d9b32SPeter Brune 
1244421d9b32SPeter Brune PetscErrorCode SNESSolve_FAS(SNES snes)
1245421d9b32SPeter Brune {
1246fa9694d7SPeter Brune   PetscErrorCode ierr;
1247fe6f9142SPeter Brune   PetscInt       i, maxits;
1248ddb5aff1SPeter Brune   Vec            X, F;
1249fe6f9142SPeter Brune   PetscReal      fnorm;
125007144faaSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
1251421d9b32SPeter Brune   PetscFunctionBegin;
1252fe6f9142SPeter Brune   maxits = snes->max_its;            /* maximum number of iterations */
1253fe6f9142SPeter Brune   snes->reason = SNES_CONVERGED_ITERATING;
1254fa9694d7SPeter Brune   X = snes->vec_sol;
1255f5a6d4f9SBarry Smith   F = snes->vec_func;
1256293a7e31SPeter Brune 
1257293a7e31SPeter Brune   /*norm setup */
1258fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1259fe6f9142SPeter Brune   snes->iter = 0;
1260fe6f9142SPeter Brune   snes->norm = 0.;
1261fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1262fe6f9142SPeter Brune   ierr = SNESComputeFunction(snes,X,F);CHKERRQ(ierr);
1263fe6f9142SPeter Brune   if (snes->domainerror) {
1264fe6f9142SPeter Brune     snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
1265fe6f9142SPeter Brune     PetscFunctionReturn(0);
1266fe6f9142SPeter Brune   }
1267fe6f9142SPeter Brune   ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1268fe6f9142SPeter Brune   if (PetscIsInfOrNanReal(fnorm)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_FP,"Infinite or not-a-number generated in norm");
1269fe6f9142SPeter Brune   ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1270fe6f9142SPeter Brune   snes->norm = fnorm;
1271fe6f9142SPeter Brune   ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1272fe6f9142SPeter Brune   SNESLogConvHistory(snes,fnorm,0);
1273fe6f9142SPeter Brune   ierr = SNESMonitor(snes,0,fnorm);CHKERRQ(ierr);
1274fe6f9142SPeter Brune 
1275fe6f9142SPeter Brune   /* set parameter for default relative tolerance convergence test */
1276fe6f9142SPeter Brune   snes->ttol = fnorm*snes->rtol;
1277fe6f9142SPeter Brune   /* test convergence */
1278fe6f9142SPeter Brune   ierr = (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1279fe6f9142SPeter Brune   if (snes->reason) PetscFunctionReturn(0);
1280fe6f9142SPeter Brune   for (i = 0; i < maxits; i++) {
1281fe6f9142SPeter Brune     /* Call general purpose update function */
1282646217ecSPeter Brune 
1283fe6f9142SPeter Brune     if (snes->ops->update) {
1284fe6f9142SPeter Brune       ierr = (*snes->ops->update)(snes, snes->iter);CHKERRQ(ierr);
1285fe6f9142SPeter Brune     }
128607144faaSPeter Brune     if (fas->fastype == SNES_FAS_MULTIPLICATIVE) {
128739bd7f45SPeter Brune       ierr = FASCycle_Multiplicative(snes, X);CHKERRQ(ierr);
128807144faaSPeter Brune     } else {
128907144faaSPeter Brune       ierr = FASCycle_Additive(snes, X);CHKERRQ(ierr);
129007144faaSPeter Brune     }
1291742fe5e2SPeter Brune 
1292742fe5e2SPeter Brune     /* check for FAS cycle divergence */
1293742fe5e2SPeter Brune     if (snes->reason != SNES_CONVERGED_ITERATING) {
1294742fe5e2SPeter Brune       PetscFunctionReturn(0);
1295742fe5e2SPeter Brune     }
1296c90fad12SPeter Brune     ierr = VecNorm(F, NORM_2, &fnorm);CHKERRQ(ierr); /* fnorm <- ||F||  */
1297c90fad12SPeter Brune     /* Monitor convergence */
1298c90fad12SPeter Brune     ierr = PetscObjectTakeAccess(snes);CHKERRQ(ierr);
1299c90fad12SPeter Brune     snes->iter = i+1;
1300c90fad12SPeter Brune     snes->norm = fnorm;
1301c90fad12SPeter Brune     ierr = PetscObjectGrantAccess(snes);CHKERRQ(ierr);
1302c90fad12SPeter Brune     SNESLogConvHistory(snes,snes->norm,0);
1303c90fad12SPeter Brune     ierr = SNESMonitor(snes,snes->iter,snes->norm);CHKERRQ(ierr);
1304c90fad12SPeter Brune     /* Test for convergence */
1305c90fad12SPeter Brune     ierr = (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);CHKERRQ(ierr);
1306c90fad12SPeter Brune     if (snes->reason) break;
1307fe6f9142SPeter Brune   }
1308fe6f9142SPeter Brune   if (i == maxits) {
1309fe6f9142SPeter Brune     ierr = PetscInfo1(snes, "Maximum number of iterations has been reached: %D\n", maxits);CHKERRQ(ierr);
1310fe6f9142SPeter Brune     if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
1311fe6f9142SPeter Brune   }
1312421d9b32SPeter Brune   PetscFunctionReturn(0);
1313421d9b32SPeter Brune }
1314