xref: /petsc/src/snes/impls/fas/fasfunc.c (revision ab8d36c96b69fea5736263697ddaa972e87bcd95)
1*ab8d36c9SPeter Brune #include "../src/snes/impls/fas/fasimpls.h" /*I  "petscsnesfas.h"  I*/
2*ab8d36c9SPeter Brune 
3*ab8d36c9SPeter Brune 
4*ab8d36c9SPeter Brune extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES *);
5*ab8d36c9SPeter Brune 
6*ab8d36c9SPeter Brune /* -------------- functions called on the fine level -------------- */
7*ab8d36c9SPeter Brune 
8*ab8d36c9SPeter Brune #undef __FUNCT__
9*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetType"
10*ab8d36c9SPeter Brune /*@
11*ab8d36c9SPeter Brune SNESFASSetType - Sets the update and correction type used for FAS.
12*ab8d36c9SPeter Brune 
13*ab8d36c9SPeter Brune Logically Collective
14*ab8d36c9SPeter Brune 
15*ab8d36c9SPeter Brune Input Parameters:
16*ab8d36c9SPeter Brune snes  - FAS context
17*ab8d36c9SPeter Brune fastype  - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE
18*ab8d36c9SPeter Brune 
19*ab8d36c9SPeter Brune .seealso: PCMGSetType()
20*ab8d36c9SPeter Brune @*/
21*ab8d36c9SPeter Brune PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
22*ab8d36c9SPeter Brune {
23*ab8d36c9SPeter Brune   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
24*ab8d36c9SPeter Brune   PetscErrorCode             ierr;
25*ab8d36c9SPeter Brune   PetscFunctionBegin;
26*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
27*ab8d36c9SPeter Brune   PetscValidLogicalCollectiveEnum(snes,fastype,2);
28*ab8d36c9SPeter Brune   fas->fastype = fastype;
29*ab8d36c9SPeter Brune   if (fas->next) ierr = SNESFASSetType(fas->next, fastype);CHKERRQ(ierr);
30*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
31*ab8d36c9SPeter Brune }
32*ab8d36c9SPeter Brune 
33*ab8d36c9SPeter Brune 
34*ab8d36c9SPeter Brune #undef __FUNCT__
35*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetType"
36*ab8d36c9SPeter Brune /*@
37*ab8d36c9SPeter Brune SNESFASGetType - Sets the update and correction type used for FAS.
38*ab8d36c9SPeter Brune 
39*ab8d36c9SPeter Brune Logically Collective
40*ab8d36c9SPeter Brune 
41*ab8d36c9SPeter Brune Input Parameters:
42*ab8d36c9SPeter Brune .snes - FAS context
43*ab8d36c9SPeter Brune 
44*ab8d36c9SPeter Brune Output Parameters:
45*ab8d36c9SPeter Brune .fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE
46*ab8d36c9SPeter Brune 
47*ab8d36c9SPeter Brune .seealso: PCMGSetType()
48*ab8d36c9SPeter Brune @*/
49*ab8d36c9SPeter Brune PetscErrorCode  SNESFASGetType(SNES snes,SNESFASType *fastype)
50*ab8d36c9SPeter Brune {
51*ab8d36c9SPeter Brune   SNES_FAS                   *fas = (SNES_FAS*)snes->data;
52*ab8d36c9SPeter Brune 
53*ab8d36c9SPeter Brune   PetscFunctionBegin;
54*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
55*ab8d36c9SPeter Brune   PetscValidPointer(fastype, 2);
56*ab8d36c9SPeter Brune   *fastype = fas->fastype;
57*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
58*ab8d36c9SPeter Brune }
59*ab8d36c9SPeter Brune 
60*ab8d36c9SPeter Brune #undef __FUNCT__
61*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetLevels"
62*ab8d36c9SPeter Brune /*@C
63*ab8d36c9SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
64*ab8d36c9SPeter Brune    Must be called before any other FAS routine.
65*ab8d36c9SPeter Brune 
66*ab8d36c9SPeter Brune    Input Parameters:
67*ab8d36c9SPeter Brune +  snes   - the snes context
68*ab8d36c9SPeter Brune .  levels - the number of levels
69*ab8d36c9SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
70*ab8d36c9SPeter Brune             problems on smaller sets of processors. Use PETSC_NULL_OBJECT for default in
71*ab8d36c9SPeter Brune             Fortran.
72*ab8d36c9SPeter Brune 
73*ab8d36c9SPeter Brune    Level: intermediate
74*ab8d36c9SPeter Brune 
75*ab8d36c9SPeter Brune    Notes:
76*ab8d36c9SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
77*ab8d36c9SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
78*ab8d36c9SPeter Brune 
79*ab8d36c9SPeter Brune .keywords: FAS, MG, set, levels, multigrid
80*ab8d36c9SPeter Brune 
81*ab8d36c9SPeter Brune .seealso: SNESFASGetLevels()
82*ab8d36c9SPeter Brune @*/
83*ab8d36c9SPeter Brune PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms) {
84*ab8d36c9SPeter Brune   PetscErrorCode ierr;
85*ab8d36c9SPeter Brune   PetscInt i;
86*ab8d36c9SPeter Brune   const char     *optionsprefix;
87*ab8d36c9SPeter Brune   char           tprefix[128];
88*ab8d36c9SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
89*ab8d36c9SPeter Brune   SNES prevsnes;
90*ab8d36c9SPeter Brune   MPI_Comm comm;
91*ab8d36c9SPeter Brune   PetscFunctionBegin;
92*ab8d36c9SPeter Brune   comm = ((PetscObject)snes)->comm;
93*ab8d36c9SPeter Brune   if (levels == fas->levels) {
94*ab8d36c9SPeter Brune     if (!comms)
95*ab8d36c9SPeter Brune       PetscFunctionReturn(0);
96*ab8d36c9SPeter Brune   }
97*ab8d36c9SPeter Brune   /* user has changed the number of levels; reset */
98*ab8d36c9SPeter Brune   ierr = SNESReset(snes);CHKERRQ(ierr);
99*ab8d36c9SPeter Brune   /* destroy any coarser levels if necessary */
100*ab8d36c9SPeter Brune   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
101*ab8d36c9SPeter Brune   fas->next = PETSC_NULL;
102*ab8d36c9SPeter Brune   fas->previous = PETSC_NULL;
103*ab8d36c9SPeter Brune   prevsnes = snes;
104*ab8d36c9SPeter Brune   /* setup the finest level */
105*ab8d36c9SPeter Brune   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
106*ab8d36c9SPeter Brune   for (i = levels - 1; i >= 0; i--) {
107*ab8d36c9SPeter Brune     if (comms) comm = comms[i];
108*ab8d36c9SPeter Brune     fas->level = i;
109*ab8d36c9SPeter Brune     fas->levels = levels;
110*ab8d36c9SPeter Brune     fas->fine = snes;
111*ab8d36c9SPeter Brune     fas->next = PETSC_NULL;
112*ab8d36c9SPeter Brune     if (i > 0) {
113*ab8d36c9SPeter Brune       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
114*ab8d36c9SPeter Brune       sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level);
115*ab8d36c9SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->next,tprefix);CHKERRQ(ierr);
116*ab8d36c9SPeter Brune       ierr = SNESAppendOptionsPrefix(fas->next,optionsprefix);CHKERRQ(ierr);
117*ab8d36c9SPeter Brune       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
118*ab8d36c9SPeter Brune       ierr = SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);CHKERRQ(ierr);
119*ab8d36c9SPeter Brune       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
120*ab8d36c9SPeter Brune       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
121*ab8d36c9SPeter Brune       prevsnes = fas->next;
122*ab8d36c9SPeter Brune       fas = (SNES_FAS *)prevsnes->data;
123*ab8d36c9SPeter Brune     }
124*ab8d36c9SPeter Brune   }
125*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
126*ab8d36c9SPeter Brune }
127*ab8d36c9SPeter Brune 
128*ab8d36c9SPeter Brune 
129*ab8d36c9SPeter Brune #undef __FUNCT__
130*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetLevels"
131*ab8d36c9SPeter Brune /*@
132*ab8d36c9SPeter Brune    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
133*ab8d36c9SPeter Brune 
134*ab8d36c9SPeter Brune    Input Parameter:
135*ab8d36c9SPeter Brune .  snes - the nonlinear solver context
136*ab8d36c9SPeter Brune 
137*ab8d36c9SPeter Brune    Output parameter:
138*ab8d36c9SPeter Brune .  levels - the number of levels
139*ab8d36c9SPeter Brune 
140*ab8d36c9SPeter Brune    Level: advanced
141*ab8d36c9SPeter Brune 
142*ab8d36c9SPeter Brune .keywords: MG, get, levels, multigrid
143*ab8d36c9SPeter Brune 
144*ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels()
145*ab8d36c9SPeter Brune @*/
146*ab8d36c9SPeter Brune PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt * levels) {
147*ab8d36c9SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
148*ab8d36c9SPeter Brune   PetscFunctionBegin;
149*ab8d36c9SPeter Brune   *levels = fas->levels;
150*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
151*ab8d36c9SPeter Brune }
152*ab8d36c9SPeter Brune 
153*ab8d36c9SPeter Brune 
154*ab8d36c9SPeter Brune #undef __FUNCT__
155*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetCycleSNES"
156*ab8d36c9SPeter Brune /*@
157*ab8d36c9SPeter Brune    SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
158*ab8d36c9SPeter Brune    level of the FAS hierarchy.
159*ab8d36c9SPeter Brune 
160*ab8d36c9SPeter Brune    Input Parameters:
161*ab8d36c9SPeter Brune +  snes    - the multigrid context
162*ab8d36c9SPeter Brune    level   - the level to get
163*ab8d36c9SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
164*ab8d36c9SPeter Brune 
165*ab8d36c9SPeter Brune    Level: advanced
166*ab8d36c9SPeter Brune 
167*ab8d36c9SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
168*ab8d36c9SPeter Brune 
169*ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
170*ab8d36c9SPeter Brune @*/
171*ab8d36c9SPeter Brune PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes) {
172*ab8d36c9SPeter Brune   SNES_FAS *fas = (SNES_FAS*)snes->data;
173*ab8d36c9SPeter Brune   PetscInt i;
174*ab8d36c9SPeter Brune 
175*ab8d36c9SPeter Brune   PetscFunctionBegin;
176*ab8d36c9SPeter Brune   if (level > fas->levels-1) SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
177*ab8d36c9SPeter Brune   if (fas->level !=  fas->levels - 1)
178*ab8d36c9SPeter Brune     SETERRQ2(((PetscObject)snes)->comm,PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level);
179*ab8d36c9SPeter Brune 
180*ab8d36c9SPeter Brune   *lsnes = snes;
181*ab8d36c9SPeter Brune   for (i = fas->level; i > level; i--) {
182*ab8d36c9SPeter Brune     *lsnes = fas->next;
183*ab8d36c9SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
184*ab8d36c9SPeter Brune   }
185*ab8d36c9SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
186*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
187*ab8d36c9SPeter Brune }
188*ab8d36c9SPeter Brune 
189*ab8d36c9SPeter Brune #undef __FUNCT__
190*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp"
191*ab8d36c9SPeter Brune /*@
192*ab8d36c9SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
193*ab8d36c9SPeter Brune    use on all levels.
194*ab8d36c9SPeter Brune 
195*ab8d36c9SPeter Brune    Logically Collective on SNES
196*ab8d36c9SPeter Brune 
197*ab8d36c9SPeter Brune    Input Parameters:
198*ab8d36c9SPeter Brune +  snes - the multigrid context
199*ab8d36c9SPeter Brune -  n    - the number of smoothing steps
200*ab8d36c9SPeter Brune 
201*ab8d36c9SPeter Brune    Options Database Key:
202*ab8d36c9SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
203*ab8d36c9SPeter Brune 
204*ab8d36c9SPeter Brune    Level: advanced
205*ab8d36c9SPeter Brune 
206*ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
207*ab8d36c9SPeter Brune 
208*ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
209*ab8d36c9SPeter Brune @*/
210*ab8d36c9SPeter Brune PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
211*ab8d36c9SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
212*ab8d36c9SPeter Brune   PetscErrorCode ierr = 0;
213*ab8d36c9SPeter Brune   PetscFunctionBegin;
214*ab8d36c9SPeter Brune   fas->max_up_it = n;
215*ab8d36c9SPeter Brune   if (!fas->smoothu) {
216*ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
217*ab8d36c9SPeter Brune   }
218*ab8d36c9SPeter Brune   ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr);
219*ab8d36c9SPeter Brune   if (fas->next) {
220*ab8d36c9SPeter Brune     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
221*ab8d36c9SPeter Brune   }
222*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
223*ab8d36c9SPeter Brune }
224*ab8d36c9SPeter Brune 
225*ab8d36c9SPeter Brune #undef __FUNCT__
226*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown"
227*ab8d36c9SPeter Brune /*@
228*ab8d36c9SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
229*ab8d36c9SPeter Brune    use on all levels.
230*ab8d36c9SPeter Brune 
231*ab8d36c9SPeter Brune    Logically Collective on SNES
232*ab8d36c9SPeter Brune 
233*ab8d36c9SPeter Brune    Input Parameters:
234*ab8d36c9SPeter Brune +  snes - the multigrid context
235*ab8d36c9SPeter Brune -  n    - the number of smoothing steps
236*ab8d36c9SPeter Brune 
237*ab8d36c9SPeter Brune    Options Database Key:
238*ab8d36c9SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
239*ab8d36c9SPeter Brune 
240*ab8d36c9SPeter Brune    Level: advanced
241*ab8d36c9SPeter Brune 
242*ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
243*ab8d36c9SPeter Brune 
244*ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
245*ab8d36c9SPeter Brune @*/
246*ab8d36c9SPeter Brune PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
247*ab8d36c9SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
248*ab8d36c9SPeter Brune   PetscErrorCode ierr = 0;
249*ab8d36c9SPeter Brune   PetscFunctionBegin;
250*ab8d36c9SPeter Brune   if (!fas->smoothu) {
251*ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);
252*ab8d36c9SPeter Brune   }
253*ab8d36c9SPeter Brune   if (!fas->smoothd) {
254*ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);
255*ab8d36c9SPeter Brune   }
256*ab8d36c9SPeter Brune   ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr);
257*ab8d36c9SPeter Brune   fas->max_down_it = n;
258*ab8d36c9SPeter Brune   if (fas->next) {
259*ab8d36c9SPeter Brune     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
260*ab8d36c9SPeter Brune   }
261*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
262*ab8d36c9SPeter Brune }
263*ab8d36c9SPeter Brune 
264*ab8d36c9SPeter Brune 
265*ab8d36c9SPeter Brune #undef __FUNCT__
266*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetCycles"
267*ab8d36c9SPeter Brune /*@
268*ab8d36c9SPeter Brune    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
269*ab8d36c9SPeter Brune    complicated cycling.
270*ab8d36c9SPeter Brune 
271*ab8d36c9SPeter Brune    Logically Collective on SNES
272*ab8d36c9SPeter Brune 
273*ab8d36c9SPeter Brune    Input Parameters:
274*ab8d36c9SPeter Brune +  snes   - the multigrid context
275*ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
276*ab8d36c9SPeter Brune 
277*ab8d36c9SPeter Brune    Options Database Key:
278*ab8d36c9SPeter Brune $  -snes_fas_cycles 1 or 2
279*ab8d36c9SPeter Brune 
280*ab8d36c9SPeter Brune    Level: advanced
281*ab8d36c9SPeter Brune 
282*ab8d36c9SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
283*ab8d36c9SPeter Brune 
284*ab8d36c9SPeter Brune .seealso: SNESFASSetCyclesOnLevel()
285*ab8d36c9SPeter Brune @*/
286*ab8d36c9SPeter Brune PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
287*ab8d36c9SPeter Brune   SNES_FAS * fas = (SNES_FAS *)snes->data;
288*ab8d36c9SPeter Brune   PetscErrorCode ierr;
289*ab8d36c9SPeter Brune   PetscBool isFine;
290*ab8d36c9SPeter Brune   PetscFunctionBegin;
291*ab8d36c9SPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);
292*ab8d36c9SPeter Brune   fas->n_cycles = cycles;
293*ab8d36c9SPeter Brune   if (!isFine)
294*ab8d36c9SPeter Brune     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
295*ab8d36c9SPeter Brune   if (fas->next) {
296*ab8d36c9SPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
297*ab8d36c9SPeter Brune   }
298*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
299*ab8d36c9SPeter Brune }
300*ab8d36c9SPeter Brune 
301*ab8d36c9SPeter Brune #undef __FUNCT__
302*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleCreateSmoother_Private"
303*ab8d36c9SPeter Brune /*
304*ab8d36c9SPeter Brune Creates the default smoother type.
305*ab8d36c9SPeter Brune 
306*ab8d36c9SPeter Brune This is SNESNRICHARDSON on each fine level and SNESLS on the coarse level.
307*ab8d36c9SPeter Brune 
308*ab8d36c9SPeter Brune  */
309*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) {
310*ab8d36c9SPeter Brune   SNES_FAS       *fas;
311*ab8d36c9SPeter Brune   const char     *optionsprefix;
312*ab8d36c9SPeter Brune   char           tprefix[128];
313*ab8d36c9SPeter Brune   PetscErrorCode ierr;
314*ab8d36c9SPeter Brune   SNES           nsmooth;
315*ab8d36c9SPeter Brune   PetscFunctionBegin;
316*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
317*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
318*ab8d36c9SPeter Brune   ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr);
319*ab8d36c9SPeter Brune   /* create the default smoother */
320*ab8d36c9SPeter Brune   ierr = SNESCreate(((PetscObject)snes)->comm, &nsmooth);CHKERRQ(ierr);
321*ab8d36c9SPeter Brune   if (fas->level == 0) {
322*ab8d36c9SPeter Brune     sprintf(tprefix,"fas_coarse_");
323*ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
324*ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
325*ab8d36c9SPeter Brune     ierr = SNESSetType(nsmooth, SNESLS);CHKERRQ(ierr);
326*ab8d36c9SPeter Brune     ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, fas->max_up_it, 1000);CHKERRQ(ierr);
327*ab8d36c9SPeter Brune     ierr = SNESSetFromOptions(nsmooth);CHKERRQ(ierr);
328*ab8d36c9SPeter Brune   } else {
329*ab8d36c9SPeter Brune     sprintf(tprefix,"fas_levels_%d_",(int)fas->level);
330*ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
331*ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
332*ab8d36c9SPeter Brune     ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr);
333*ab8d36c9SPeter Brune     ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, fas->max_up_it, 1000);CHKERRQ(ierr);
334*ab8d36c9SPeter Brune     ierr = SNESSetFromOptions(nsmooth);CHKERRQ(ierr);
335*ab8d36c9SPeter Brune   }
336*ab8d36c9SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
337*ab8d36c9SPeter Brune   *smooth = nsmooth;
338*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
339*ab8d36c9SPeter Brune }
340*ab8d36c9SPeter Brune 
341*ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */
342*ab8d36c9SPeter Brune 
343*ab8d36c9SPeter Brune #undef __FUNCT__
344*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleSetCycles"
345*ab8d36c9SPeter Brune /*@
346*ab8d36c9SPeter Brune    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
347*ab8d36c9SPeter Brune 
348*ab8d36c9SPeter Brune    Logically Collective on SNES
349*ab8d36c9SPeter Brune 
350*ab8d36c9SPeter Brune    Input Parameters:
351*ab8d36c9SPeter Brune +  snes   - the multigrid context
352*ab8d36c9SPeter Brune .  level  - the level to set the number of cycles on
353*ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
354*ab8d36c9SPeter Brune 
355*ab8d36c9SPeter Brune    Level: advanced
356*ab8d36c9SPeter Brune 
357*ab8d36c9SPeter Brune .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid
358*ab8d36c9SPeter Brune 
359*ab8d36c9SPeter Brune .seealso: SNESFASSetCycles()
360*ab8d36c9SPeter Brune @*/
361*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) {
362*ab8d36c9SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
363*ab8d36c9SPeter Brune   PetscErrorCode ierr;
364*ab8d36c9SPeter Brune   PetscFunctionBegin;
365*ab8d36c9SPeter Brune   fas->n_cycles = cycles;
366*ab8d36c9SPeter Brune   ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
367*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
368*ab8d36c9SPeter Brune }
369*ab8d36c9SPeter Brune 
370*ab8d36c9SPeter Brune 
371*ab8d36c9SPeter Brune #undef __FUNCT__
372*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmoother"
373*ab8d36c9SPeter Brune /*@
374*ab8d36c9SPeter Brune    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
375*ab8d36c9SPeter Brune 
376*ab8d36c9SPeter Brune    Logically Collective on SNES
377*ab8d36c9SPeter Brune 
378*ab8d36c9SPeter Brune    Input Parameters:
379*ab8d36c9SPeter Brune .  snes   - the multigrid context
380*ab8d36c9SPeter Brune 
381*ab8d36c9SPeter Brune    Output Parameters:
382*ab8d36c9SPeter Brune .  smooth - the smoother
383*ab8d36c9SPeter Brune 
384*ab8d36c9SPeter Brune    Level: advanced
385*ab8d36c9SPeter Brune 
386*ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
387*ab8d36c9SPeter Brune 
388*ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
389*ab8d36c9SPeter Brune @*/
390*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
391*ab8d36c9SPeter Brune {
392*ab8d36c9SPeter Brune   SNES_FAS       *fas;
393*ab8d36c9SPeter Brune   PetscFunctionBegin;
394*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
395*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
396*ab8d36c9SPeter Brune   *smooth = fas->smoothd;
397*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
398*ab8d36c9SPeter Brune }
399*ab8d36c9SPeter Brune #undef __FUNCT__
400*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherUp"
401*ab8d36c9SPeter Brune /*@
402*ab8d36c9SPeter Brune    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
403*ab8d36c9SPeter Brune 
404*ab8d36c9SPeter Brune    Logically Collective on SNES
405*ab8d36c9SPeter Brune 
406*ab8d36c9SPeter Brune    Input Parameters:
407*ab8d36c9SPeter Brune .  snes   - the multigrid context
408*ab8d36c9SPeter Brune 
409*ab8d36c9SPeter Brune    Output Parameters:
410*ab8d36c9SPeter Brune .  smoothu - the smoother
411*ab8d36c9SPeter Brune 
412*ab8d36c9SPeter Brune    Notes:
413*ab8d36c9SPeter Brune    Returns the downsmoother if no up smoother is available.  This enables transparent
414*ab8d36c9SPeter Brune    default behavior in the process of the solve.
415*ab8d36c9SPeter Brune 
416*ab8d36c9SPeter Brune    Level: advanced
417*ab8d36c9SPeter Brune 
418*ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
419*ab8d36c9SPeter Brune 
420*ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
421*ab8d36c9SPeter Brune @*/
422*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
423*ab8d36c9SPeter Brune {
424*ab8d36c9SPeter Brune   SNES_FAS       *fas;
425*ab8d36c9SPeter Brune   PetscFunctionBegin;
426*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
427*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
428*ab8d36c9SPeter Brune   if (!fas->smoothu) {
429*ab8d36c9SPeter Brune     *smoothu = fas->smoothd;
430*ab8d36c9SPeter Brune   } else {
431*ab8d36c9SPeter Brune     *smoothu = fas->smoothu;
432*ab8d36c9SPeter Brune   }
433*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
434*ab8d36c9SPeter Brune }
435*ab8d36c9SPeter Brune 
436*ab8d36c9SPeter Brune #undef __FUNCT__
437*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherDown"
438*ab8d36c9SPeter Brune /*@
439*ab8d36c9SPeter Brune    SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
440*ab8d36c9SPeter Brune 
441*ab8d36c9SPeter Brune    Logically Collective on SNES
442*ab8d36c9SPeter Brune 
443*ab8d36c9SPeter Brune    Input Parameters:
444*ab8d36c9SPeter Brune .  snes   - the multigrid context
445*ab8d36c9SPeter Brune 
446*ab8d36c9SPeter Brune    Output Parameters:
447*ab8d36c9SPeter Brune .  smoothd - the smoother
448*ab8d36c9SPeter Brune 
449*ab8d36c9SPeter Brune    Level: advanced
450*ab8d36c9SPeter Brune 
451*ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
452*ab8d36c9SPeter Brune 
453*ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
454*ab8d36c9SPeter Brune @*/
455*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
456*ab8d36c9SPeter Brune {
457*ab8d36c9SPeter Brune   SNES_FAS       *fas;
458*ab8d36c9SPeter Brune   PetscFunctionBegin;
459*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
460*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
461*ab8d36c9SPeter Brune   *smoothd = fas->smoothd;
462*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
463*ab8d36c9SPeter Brune }
464*ab8d36c9SPeter Brune 
465*ab8d36c9SPeter Brune 
466*ab8d36c9SPeter Brune #undef __FUNCT__
467*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetCorrection"
468*ab8d36c9SPeter Brune /*@
469*ab8d36c9SPeter Brune    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level
470*ab8d36c9SPeter Brune 
471*ab8d36c9SPeter Brune    Logically Collective on SNES
472*ab8d36c9SPeter Brune 
473*ab8d36c9SPeter Brune    Input Parameters:
474*ab8d36c9SPeter Brune .  snes   - the multigrid context
475*ab8d36c9SPeter Brune 
476*ab8d36c9SPeter Brune    Output Parameters:
477*ab8d36c9SPeter Brune .  correction - the coarse correction on this level
478*ab8d36c9SPeter Brune 
479*ab8d36c9SPeter Brune    Notes:
480*ab8d36c9SPeter Brune    Returns PETSC_NULL on the coarsest level.
481*ab8d36c9SPeter Brune 
482*ab8d36c9SPeter Brune    Level: advanced
483*ab8d36c9SPeter Brune 
484*ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
485*ab8d36c9SPeter Brune 
486*ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
487*ab8d36c9SPeter Brune @*/
488*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
489*ab8d36c9SPeter Brune {
490*ab8d36c9SPeter Brune   SNES_FAS       *fas;
491*ab8d36c9SPeter Brune   PetscFunctionBegin;
492*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
493*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
494*ab8d36c9SPeter Brune   *correction = fas->next;
495*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
496*ab8d36c9SPeter Brune }
497*ab8d36c9SPeter Brune 
498*ab8d36c9SPeter Brune #undef __FUNCT__
499*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInterpolation"
500*ab8d36c9SPeter Brune /*@
501*ab8d36c9SPeter Brune    SNESFASCycleGetInterpolation - Gets the interpolation on this level
502*ab8d36c9SPeter Brune 
503*ab8d36c9SPeter Brune    Logically Collective on SNES
504*ab8d36c9SPeter Brune 
505*ab8d36c9SPeter Brune    Input Parameters:
506*ab8d36c9SPeter Brune .  snes   - the multigrid context
507*ab8d36c9SPeter Brune 
508*ab8d36c9SPeter Brune    Output Parameters:
509*ab8d36c9SPeter Brune .  mat    - the interpolation operator on this level
510*ab8d36c9SPeter Brune 
511*ab8d36c9SPeter Brune    Level: developer
512*ab8d36c9SPeter Brune 
513*ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
514*ab8d36c9SPeter Brune 
515*ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
516*ab8d36c9SPeter Brune @*/
517*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
518*ab8d36c9SPeter Brune {
519*ab8d36c9SPeter Brune   SNES_FAS       *fas;
520*ab8d36c9SPeter Brune   PetscFunctionBegin;
521*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
522*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
523*ab8d36c9SPeter Brune   *mat = fas->interpolate;
524*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
525*ab8d36c9SPeter Brune }
526*ab8d36c9SPeter Brune 
527*ab8d36c9SPeter Brune 
528*ab8d36c9SPeter Brune #undef __FUNCT__
529*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRestriction"
530*ab8d36c9SPeter Brune /*@
531*ab8d36c9SPeter Brune    SNESFASCycleGetRestriction - Gets the restriction on this level
532*ab8d36c9SPeter Brune 
533*ab8d36c9SPeter Brune    Logically Collective on SNES
534*ab8d36c9SPeter Brune 
535*ab8d36c9SPeter Brune    Input Parameters:
536*ab8d36c9SPeter Brune .  snes   - the multigrid context
537*ab8d36c9SPeter Brune 
538*ab8d36c9SPeter Brune    Output Parameters:
539*ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
540*ab8d36c9SPeter Brune 
541*ab8d36c9SPeter Brune    Level: developer
542*ab8d36c9SPeter Brune 
543*ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
544*ab8d36c9SPeter Brune 
545*ab8d36c9SPeter Brune .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
546*ab8d36c9SPeter Brune @*/
547*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
548*ab8d36c9SPeter Brune {
549*ab8d36c9SPeter Brune   SNES_FAS       *fas;
550*ab8d36c9SPeter Brune   PetscFunctionBegin;
551*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
552*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
553*ab8d36c9SPeter Brune   *mat = fas->restrct;
554*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
555*ab8d36c9SPeter Brune }
556*ab8d36c9SPeter Brune 
557*ab8d36c9SPeter Brune 
558*ab8d36c9SPeter Brune #undef __FUNCT__
559*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInjection"
560*ab8d36c9SPeter Brune /*@
561*ab8d36c9SPeter Brune    SNESFASCycleGetInjection - Gets the injection on this level
562*ab8d36c9SPeter Brune 
563*ab8d36c9SPeter Brune    Logically Collective on SNES
564*ab8d36c9SPeter Brune 
565*ab8d36c9SPeter Brune    Input Parameters:
566*ab8d36c9SPeter Brune .  snes   - the multigrid context
567*ab8d36c9SPeter Brune 
568*ab8d36c9SPeter Brune    Output Parameters:
569*ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
570*ab8d36c9SPeter Brune 
571*ab8d36c9SPeter Brune    Level: developer
572*ab8d36c9SPeter Brune 
573*ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
574*ab8d36c9SPeter Brune 
575*ab8d36c9SPeter Brune .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
576*ab8d36c9SPeter Brune @*/
577*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
578*ab8d36c9SPeter Brune {
579*ab8d36c9SPeter Brune   SNES_FAS       *fas;
580*ab8d36c9SPeter Brune   PetscFunctionBegin;
581*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
582*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
583*ab8d36c9SPeter Brune   *mat = fas->inject;
584*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
585*ab8d36c9SPeter Brune }
586*ab8d36c9SPeter Brune 
587*ab8d36c9SPeter Brune #undef __FUNCT__
588*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRScale"
589*ab8d36c9SPeter Brune /*@
590*ab8d36c9SPeter Brune    SNESFASCycleGetRScale - Gets the injection on this level
591*ab8d36c9SPeter Brune 
592*ab8d36c9SPeter Brune    Logically Collective on SNES
593*ab8d36c9SPeter Brune 
594*ab8d36c9SPeter Brune    Input Parameters:
595*ab8d36c9SPeter Brune .  snes   - the multigrid context
596*ab8d36c9SPeter Brune 
597*ab8d36c9SPeter Brune    Output Parameters:
598*ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
599*ab8d36c9SPeter Brune 
600*ab8d36c9SPeter Brune    Level: developer
601*ab8d36c9SPeter Brune 
602*ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
603*ab8d36c9SPeter Brune 
604*ab8d36c9SPeter Brune .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
605*ab8d36c9SPeter Brune @*/
606*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
607*ab8d36c9SPeter Brune {
608*ab8d36c9SPeter Brune   SNES_FAS       *fas;
609*ab8d36c9SPeter Brune   PetscFunctionBegin;
610*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
611*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
612*ab8d36c9SPeter Brune   *vec = fas->rscale;
613*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
614*ab8d36c9SPeter Brune }
615*ab8d36c9SPeter Brune 
616*ab8d36c9SPeter Brune #undef __FUNCT__
617*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleIsFine"
618*ab8d36c9SPeter Brune /*@
619*ab8d36c9SPeter Brune    SNESFASCycleIsFine - Determines if a given cycle is the fine level.
620*ab8d36c9SPeter Brune 
621*ab8d36c9SPeter Brune    Logically Collective on SNES
622*ab8d36c9SPeter Brune 
623*ab8d36c9SPeter Brune    Input Parameters:
624*ab8d36c9SPeter Brune .  snes   - the FAS context
625*ab8d36c9SPeter Brune 
626*ab8d36c9SPeter Brune    Output Parameters:
627*ab8d36c9SPeter Brune .  flg - indicates if this is the fine level or not
628*ab8d36c9SPeter Brune 
629*ab8d36c9SPeter Brune    Level: advanced
630*ab8d36c9SPeter Brune 
631*ab8d36c9SPeter Brune .keywords: SNES, FAS
632*ab8d36c9SPeter Brune 
633*ab8d36c9SPeter Brune .seealso: SNESFASSetLevels()
634*ab8d36c9SPeter Brune @*/
635*ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
636*ab8d36c9SPeter Brune {
637*ab8d36c9SPeter Brune   SNES_FAS       *fas;
638*ab8d36c9SPeter Brune   PetscFunctionBegin;
639*ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
640*ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
641*ab8d36c9SPeter Brune   if (fas->level == fas->levels - 1) {
642*ab8d36c9SPeter Brune     *flg = PETSC_TRUE;
643*ab8d36c9SPeter Brune   } else {
644*ab8d36c9SPeter Brune     *flg = PETSC_FALSE;
645*ab8d36c9SPeter Brune   }
646*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
647*ab8d36c9SPeter Brune }
648*ab8d36c9SPeter Brune 
649*ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */
650*ab8d36c9SPeter Brune 
651*ab8d36c9SPeter Brune #undef __FUNCT__
652*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
653*ab8d36c9SPeter Brune /*@
654*ab8d36c9SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
655*ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
656*ab8d36c9SPeter Brune 
657*ab8d36c9SPeter Brune    Input Parameters:
658*ab8d36c9SPeter Brune +  snes      - the multigrid context
659*ab8d36c9SPeter Brune .  mat       - the interpolation operator
660*ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
661*ab8d36c9SPeter Brune 
662*ab8d36c9SPeter Brune    Level: advanced
663*ab8d36c9SPeter Brune 
664*ab8d36c9SPeter Brune    Notes:
665*ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the restriction
666*ab8d36c9SPeter Brune     for the same level.
667*ab8d36c9SPeter Brune 
668*ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
669*ab8d36c9SPeter Brune     out from the matrix size which one it is.
670*ab8d36c9SPeter Brune 
671*ab8d36c9SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
672*ab8d36c9SPeter Brune 
673*ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
674*ab8d36c9SPeter Brune @*/
675*ab8d36c9SPeter Brune PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
676*ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
677*ab8d36c9SPeter Brune   PetscErrorCode ierr;
678*ab8d36c9SPeter Brune   SNES           levelsnes;
679*ab8d36c9SPeter Brune   PetscFunctionBegin;
680*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
681*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
682*ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
683*ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
684*ab8d36c9SPeter Brune   fas->interpolate = mat;
685*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
686*ab8d36c9SPeter Brune }
687*ab8d36c9SPeter Brune 
688*ab8d36c9SPeter Brune #undef __FUNCT__
689*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInterpolation"
690*ab8d36c9SPeter Brune /*@
691*ab8d36c9SPeter Brune    SNESFASGetInterpolation - Gets the matrix used to calculate the
692*ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
693*ab8d36c9SPeter Brune 
694*ab8d36c9SPeter Brune    Input Parameters:
695*ab8d36c9SPeter Brune +  snes      - the multigrid context
696*ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
697*ab8d36c9SPeter Brune 
698*ab8d36c9SPeter Brune    Output Parameters:
699*ab8d36c9SPeter Brune .  mat       - the interpolation operator
700*ab8d36c9SPeter Brune 
701*ab8d36c9SPeter Brune    Level: advanced
702*ab8d36c9SPeter Brune 
703*ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, interpolate, level
704*ab8d36c9SPeter Brune 
705*ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
706*ab8d36c9SPeter Brune @*/
707*ab8d36c9SPeter Brune PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) {
708*ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
709*ab8d36c9SPeter Brune   PetscErrorCode ierr;
710*ab8d36c9SPeter Brune   SNES           levelsnes;
711*ab8d36c9SPeter Brune   PetscFunctionBegin;
712*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
713*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
714*ab8d36c9SPeter Brune   *mat = fas->interpolate;
715*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
716*ab8d36c9SPeter Brune }
717*ab8d36c9SPeter Brune 
718*ab8d36c9SPeter Brune #undef __FUNCT__
719*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
720*ab8d36c9SPeter Brune /*@
721*ab8d36c9SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
722*ab8d36c9SPeter Brune    from level l to l-1.
723*ab8d36c9SPeter Brune 
724*ab8d36c9SPeter Brune    Input Parameters:
725*ab8d36c9SPeter Brune +  snes  - the multigrid context
726*ab8d36c9SPeter Brune .  mat   - the restriction matrix
727*ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
728*ab8d36c9SPeter Brune 
729*ab8d36c9SPeter Brune    Level: advanced
730*ab8d36c9SPeter Brune 
731*ab8d36c9SPeter Brune    Notes:
732*ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the interpolation
733*ab8d36c9SPeter Brune     for the same level.
734*ab8d36c9SPeter Brune 
735*ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
736*ab8d36c9SPeter Brune     out from the matrix size which one it is.
737*ab8d36c9SPeter Brune 
738*ab8d36c9SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
739*ab8d36c9SPeter Brune     is used.
740*ab8d36c9SPeter Brune 
741*ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
742*ab8d36c9SPeter Brune 
743*ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
744*ab8d36c9SPeter Brune @*/
745*ab8d36c9SPeter Brune PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
746*ab8d36c9SPeter Brune   SNES_FAS * fas =  (SNES_FAS *)snes->data;
747*ab8d36c9SPeter Brune   PetscErrorCode ierr;
748*ab8d36c9SPeter Brune   SNES levelsnes;
749*ab8d36c9SPeter Brune   PetscFunctionBegin;
750*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
751*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
752*ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
753*ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
754*ab8d36c9SPeter Brune   fas->restrct = mat;
755*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
756*ab8d36c9SPeter Brune }
757*ab8d36c9SPeter Brune 
758*ab8d36c9SPeter Brune #undef __FUNCT__
759*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetRestriction"
760*ab8d36c9SPeter Brune /*@
761*ab8d36c9SPeter Brune    SNESFASGetRestriction - Gets the matrix used to calculate the
762*ab8d36c9SPeter Brune    restriction from l to the l-1th level
763*ab8d36c9SPeter Brune 
764*ab8d36c9SPeter Brune    Input Parameters:
765*ab8d36c9SPeter Brune +  snes      - the multigrid context
766*ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
767*ab8d36c9SPeter Brune 
768*ab8d36c9SPeter Brune    Output Parameters:
769*ab8d36c9SPeter Brune .  mat       - the interpolation operator
770*ab8d36c9SPeter Brune 
771*ab8d36c9SPeter Brune    Level: advanced
772*ab8d36c9SPeter Brune 
773*ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, restrict, level
774*ab8d36c9SPeter Brune 
775*ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
776*ab8d36c9SPeter Brune @*/
777*ab8d36c9SPeter Brune PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) {
778*ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
779*ab8d36c9SPeter Brune   PetscErrorCode ierr;
780*ab8d36c9SPeter Brune   SNES           levelsnes;
781*ab8d36c9SPeter Brune   PetscFunctionBegin;
782*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
783*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
784*ab8d36c9SPeter Brune   *mat = fas->restrct;
785*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
786*ab8d36c9SPeter Brune }
787*ab8d36c9SPeter Brune 
788*ab8d36c9SPeter Brune 
789*ab8d36c9SPeter Brune #undef __FUNCT__
790*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInjection"
791*ab8d36c9SPeter Brune /*@
792*ab8d36c9SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
793*ab8d36c9SPeter Brune    from level l to l-1.
794*ab8d36c9SPeter Brune 
795*ab8d36c9SPeter Brune    Input Parameters:
796*ab8d36c9SPeter Brune  +  snes  - the multigrid context
797*ab8d36c9SPeter Brune .  mat   - the restriction matrix
798*ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
799*ab8d36c9SPeter Brune 
800*ab8d36c9SPeter Brune    Level: advanced
801*ab8d36c9SPeter Brune 
802*ab8d36c9SPeter Brune    Notes:
803*ab8d36c9SPeter Brune          If you do not set this, the restriction and rscale is used to
804*ab8d36c9SPeter Brune    project the solution instead.
805*ab8d36c9SPeter Brune 
806*ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
807*ab8d36c9SPeter Brune 
808*ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
809*ab8d36c9SPeter Brune @*/
810*ab8d36c9SPeter Brune PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
811*ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
812*ab8d36c9SPeter Brune   PetscErrorCode ierr;
813*ab8d36c9SPeter Brune   SNES           levelsnes;
814*ab8d36c9SPeter Brune   PetscFunctionBegin;
815*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
816*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
817*ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
818*ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
819*ab8d36c9SPeter Brune   fas->inject = mat;
820*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
821*ab8d36c9SPeter Brune }
822*ab8d36c9SPeter Brune 
823*ab8d36c9SPeter Brune 
824*ab8d36c9SPeter Brune #undef __FUNCT__
825*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInjection"
826*ab8d36c9SPeter Brune /*@
827*ab8d36c9SPeter Brune    SNESFASGetInjection - Gets the matrix used to calculate the
828*ab8d36c9SPeter Brune    injection from l-1 to the lth level
829*ab8d36c9SPeter Brune 
830*ab8d36c9SPeter Brune    Input Parameters:
831*ab8d36c9SPeter Brune +  snes      - the multigrid context
832*ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
833*ab8d36c9SPeter Brune 
834*ab8d36c9SPeter Brune    Output Parameters:
835*ab8d36c9SPeter Brune .  mat       - the injection operator
836*ab8d36c9SPeter Brune 
837*ab8d36c9SPeter Brune    Level: advanced
838*ab8d36c9SPeter Brune 
839*ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, restrict, level
840*ab8d36c9SPeter Brune 
841*ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
842*ab8d36c9SPeter Brune @*/
843*ab8d36c9SPeter Brune PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) {
844*ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
845*ab8d36c9SPeter Brune   PetscErrorCode ierr;
846*ab8d36c9SPeter Brune   SNES           levelsnes;
847*ab8d36c9SPeter Brune   PetscFunctionBegin;
848*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
849*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
850*ab8d36c9SPeter Brune   *mat = fas->inject;
851*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
852*ab8d36c9SPeter Brune }
853*ab8d36c9SPeter Brune 
854*ab8d36c9SPeter Brune #undef __FUNCT__
855*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
856*ab8d36c9SPeter Brune /*@
857*ab8d36c9SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
858*ab8d36c9SPeter Brune    operator from level l to l-1.
859*ab8d36c9SPeter Brune 
860*ab8d36c9SPeter Brune    Input Parameters:
861*ab8d36c9SPeter Brune +  snes   - the multigrid context
862*ab8d36c9SPeter Brune .  rscale - the restriction scaling
863*ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
864*ab8d36c9SPeter Brune 
865*ab8d36c9SPeter Brune    Level: advanced
866*ab8d36c9SPeter Brune 
867*ab8d36c9SPeter Brune    Notes:
868*ab8d36c9SPeter Brune          This is only used in the case that the injection is not set.
869*ab8d36c9SPeter Brune 
870*ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
871*ab8d36c9SPeter Brune 
872*ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
873*ab8d36c9SPeter Brune @*/
874*ab8d36c9SPeter Brune PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
875*ab8d36c9SPeter Brune   SNES_FAS * fas;
876*ab8d36c9SPeter Brune   PetscErrorCode ierr;
877*ab8d36c9SPeter Brune   SNES levelsnes;
878*ab8d36c9SPeter Brune   PetscFunctionBegin;
879*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
880*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
881*ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
882*ab8d36c9SPeter Brune   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
883*ab8d36c9SPeter Brune   fas->rscale = rscale;
884*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
885*ab8d36c9SPeter Brune }
886*ab8d36c9SPeter Brune 
887*ab8d36c9SPeter Brune #undef __FUNCT__
888*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmoother"
889*ab8d36c9SPeter Brune /*@
890*ab8d36c9SPeter Brune    SNESFASGetSmoother - Gets the default smoother on a level.
891*ab8d36c9SPeter Brune 
892*ab8d36c9SPeter Brune    Input Parameters:
893*ab8d36c9SPeter Brune +  snes   - the multigrid context
894*ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
895*ab8d36c9SPeter Brune 
896*ab8d36c9SPeter Brune    Output Parameters:
897*ab8d36c9SPeter Brune    smooth  - the smoother
898*ab8d36c9SPeter Brune 
899*ab8d36c9SPeter Brune    Level: advanced
900*ab8d36c9SPeter Brune 
901*ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
902*ab8d36c9SPeter Brune 
903*ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
904*ab8d36c9SPeter Brune @*/
905*ab8d36c9SPeter Brune PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) {
906*ab8d36c9SPeter Brune   SNES_FAS * fas;
907*ab8d36c9SPeter Brune   PetscErrorCode ierr;
908*ab8d36c9SPeter Brune   SNES levelsnes;
909*ab8d36c9SPeter Brune   PetscFunctionBegin;
910*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
911*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
912*ab8d36c9SPeter Brune   if (!fas->smoothd) {
913*ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
914*ab8d36c9SPeter Brune   }
915*ab8d36c9SPeter Brune   *smooth = fas->smoothd;
916*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
917*ab8d36c9SPeter Brune }
918*ab8d36c9SPeter Brune 
919*ab8d36c9SPeter Brune #undef __FUNCT__
920*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherDown"
921*ab8d36c9SPeter Brune /*@
922*ab8d36c9SPeter Brune    SNESFASGetSmootherDown - Gets the downsmoother on a level.
923*ab8d36c9SPeter Brune 
924*ab8d36c9SPeter Brune    Input Parameters:
925*ab8d36c9SPeter Brune +  snes   - the multigrid context
926*ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
927*ab8d36c9SPeter Brune 
928*ab8d36c9SPeter Brune    Output Parameters:
929*ab8d36c9SPeter Brune    smooth  - the smoother
930*ab8d36c9SPeter Brune 
931*ab8d36c9SPeter Brune    Level: advanced
932*ab8d36c9SPeter Brune 
933*ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
934*ab8d36c9SPeter Brune 
935*ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
936*ab8d36c9SPeter Brune @*/
937*ab8d36c9SPeter Brune PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) {
938*ab8d36c9SPeter Brune   SNES_FAS * fas;
939*ab8d36c9SPeter Brune   PetscErrorCode ierr;
940*ab8d36c9SPeter Brune   SNES levelsnes;
941*ab8d36c9SPeter Brune   PetscFunctionBegin;
942*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
943*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
944*ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
945*ab8d36c9SPeter Brune   if (!fas->smoothd) {
946*ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
947*ab8d36c9SPeter Brune   }
948*ab8d36c9SPeter Brune   if (!fas->smoothu) {
949*ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
950*ab8d36c9SPeter Brune   }
951*ab8d36c9SPeter Brune   *smooth = fas->smoothd;
952*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
953*ab8d36c9SPeter Brune }
954*ab8d36c9SPeter Brune 
955*ab8d36c9SPeter Brune #undef __FUNCT__
956*ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherUp"
957*ab8d36c9SPeter Brune /*@
958*ab8d36c9SPeter Brune    SNESFASGetSmootherUp - Gets the upsmoother on a level.
959*ab8d36c9SPeter Brune 
960*ab8d36c9SPeter Brune    Input Parameters:
961*ab8d36c9SPeter Brune +  snes   - the multigrid context
962*ab8d36c9SPeter Brune -  level  - the level (0 is coarsest)
963*ab8d36c9SPeter Brune 
964*ab8d36c9SPeter Brune    Output Parameters:
965*ab8d36c9SPeter Brune    smooth  - the smoother
966*ab8d36c9SPeter Brune 
967*ab8d36c9SPeter Brune    Level: advanced
968*ab8d36c9SPeter Brune 
969*ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
970*ab8d36c9SPeter Brune 
971*ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
972*ab8d36c9SPeter Brune @*/
973*ab8d36c9SPeter Brune PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) {
974*ab8d36c9SPeter Brune   SNES_FAS * fas;
975*ab8d36c9SPeter Brune   PetscErrorCode ierr;
976*ab8d36c9SPeter Brune   SNES levelsnes;
977*ab8d36c9SPeter Brune   PetscFunctionBegin;
978*ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
979*ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
980*ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
981*ab8d36c9SPeter Brune   if (!fas->smoothd) {
982*ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
983*ab8d36c9SPeter Brune   }
984*ab8d36c9SPeter Brune   if (!fas->smoothu) {
985*ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
986*ab8d36c9SPeter Brune   }
987*ab8d36c9SPeter Brune   *smooth = fas->smoothu;
988*ab8d36c9SPeter Brune   PetscFunctionReturn(0);
989*ab8d36c9SPeter Brune }
990