xref: /petsc/src/snes/impls/fas/fasfunc.c (revision c094ef4021e955ef5f85f7d8a1bbc6ed64ba7621)
1 #include <../src/snes/impls/fas/fasimpls.h> /*I  "petscsnes.h"  I*/
2 
3 
4 extern PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES, SNES*);
5 
6 /* -------------- functions called on the fine level -------------- */
7 
8 #undef __FUNCT__
9 #define __FUNCT__ "SNESFASSetType"
10 /*@
11     SNESFASSetType - Sets the update and correction type used for FAS.
12 
13    Logically Collective
14 
15 Input Parameters:
16 + snes  - FAS context
17 - fastype  - SNES_FAS_ADDITIVE, SNES_FAS_MULTIPLICATIVE, SNES_FAS_FULL, or SNES_FAS_KASKADE
18 
19 Level: intermediate
20 
21 .seealso: PCMGSetType()
22 @*/
23 PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
24 {
25   SNES_FAS       *fas = (SNES_FAS*)snes->data;
26   PetscErrorCode ierr;
27 
28   PetscFunctionBegin;
29   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
30   PetscValidLogicalCollectiveEnum(snes,fastype,2);
31   fas->fastype = fastype;
32   if (fas->next) {
33     ierr = SNESFASSetType(fas->next, fastype);CHKERRQ(ierr);
34   }
35   PetscFunctionReturn(0);
36 }
37 
38 
39 #undef __FUNCT__
40 #define __FUNCT__ "SNESFASGetType"
41 /*@
42 SNESFASGetType - Sets the update and correction type used for FAS.
43 
44 Logically Collective
45 
46 Input Parameters:
47 . snes - FAS context
48 
49 Output Parameters:
50 . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE
51 
52 Level: intermediate
53 
54 .seealso: PCMGSetType()
55 @*/
56 PetscErrorCode  SNESFASGetType(SNES snes,SNESFASType *fastype)
57 {
58   SNES_FAS *fas = (SNES_FAS*)snes->data;
59 
60   PetscFunctionBegin;
61   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
62   PetscValidPointer(fastype, 2);
63   *fastype = fas->fastype;
64   PetscFunctionReturn(0);
65 }
66 
67 #undef __FUNCT__
68 #define __FUNCT__ "SNESFASSetLevels"
69 /*@C
70    SNESFASSetLevels - Sets the number of levels to use with FAS.
71    Must be called before any other FAS routine.
72 
73    Input Parameters:
74 +  snes   - the snes context
75 .  levels - the number of levels
76 -  comms  - optional communicators for each level; this is to allow solving the coarser
77             problems on smaller sets of processors. Use NULL_OBJECT for default in
78             Fortran.
79 
80    Level: intermediate
81 
82    Notes:
83      If the number of levels is one then the multigrid uses the -fas_levels prefix
84   for setting the level options rather than the -fas_coarse prefix.
85 
86 .keywords: FAS, MG, set, levels, multigrid
87 
88 .seealso: SNESFASGetLevels()
89 @*/
90 PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm * comms)
91 {
92   PetscErrorCode ierr;
93   PetscInt       i;
94   const char     *optionsprefix;
95   char           tprefix[128];
96   SNES_FAS       *fas = (SNES_FAS*)snes->data;
97   SNES           prevsnes;
98   MPI_Comm       comm;
99 
100   PetscFunctionBegin;
101   ierr = PetscObjectGetComm((PetscObject)snes,&comm);CHKERRQ(ierr);
102   if (levels == fas->levels) {
103     if (!comms) PetscFunctionReturn(0);
104   }
105   /* user has changed the number of levels; reset */
106   ierr = SNESReset(snes);CHKERRQ(ierr);
107   /* destroy any coarser levels if necessary */
108   if (fas->next) SNESDestroy(&fas->next);CHKERRQ(ierr);
109   fas->next     = NULL;
110   fas->previous = NULL;
111   prevsnes      = snes;
112   /* setup the finest level */
113   ierr = SNESGetOptionsPrefix(snes, &optionsprefix);CHKERRQ(ierr);
114   ierr = PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1);CHKERRQ(ierr);
115   for (i = levels - 1; i >= 0; i--) {
116     if (comms) comm = comms[i];
117     fas->level  = i;
118     fas->levels = levels;
119     fas->fine   = snes;
120     fas->next   = NULL;
121     if (i > 0) {
122       ierr = SNESCreate(comm, &fas->next);CHKERRQ(ierr);
123       ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr);
124       sprintf(tprefix,"fas_levels_%d_cycle_",(int)fas->level);
125       ierr = SNESAppendOptionsPrefix(fas->next,optionsprefix);CHKERRQ(ierr);
126       ierr = SNESAppendOptionsPrefix(fas->next,tprefix);CHKERRQ(ierr);
127       ierr = SNESSetType(fas->next, SNESFAS);CHKERRQ(ierr);
128       ierr = SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs);CHKERRQ(ierr);
129       ierr = PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i);CHKERRQ(ierr);
130       ierr = PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1);CHKERRQ(ierr);
131 
132       ((SNES_FAS*)fas->next->data)->previous = prevsnes;
133 
134       prevsnes = fas->next;
135       fas      = (SNES_FAS*)prevsnes->data;
136     }
137   }
138   PetscFunctionReturn(0);
139 }
140 
141 
142 #undef __FUNCT__
143 #define __FUNCT__ "SNESFASGetLevels"
144 /*@
145    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
146 
147    Input Parameter:
148 .  snes - the nonlinear solver context
149 
150    Output parameter:
151 .  levels - the number of levels
152 
153    Level: advanced
154 
155 .keywords: MG, get, levels, multigrid
156 
157 .seealso: SNESFASSetLevels(), PCMGGetLevels()
158 @*/
159 PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
160 {
161   SNES_FAS * fas = (SNES_FAS*)snes->data;
162 
163   PetscFunctionBegin;
164   *levels = fas->levels;
165   PetscFunctionReturn(0);
166 }
167 
168 
169 #undef __FUNCT__
170 #define __FUNCT__ "SNESFASGetCycleSNES"
171 /*@
172    SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
173    level of the FAS hierarchy.
174 
175    Input Parameters:
176 +  snes    - the multigrid context
177    level   - the level to get
178 -  lsnes   - whether to use the nonlinear smoother or not
179 
180    Level: advanced
181 
182 .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
183 
184 .seealso: SNESFASSetLevels(), SNESFASGetLevels()
185 @*/
186 PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes)
187 {
188   SNES_FAS *fas = (SNES_FAS*)snes->data;
189   PetscInt i;
190 
191   PetscFunctionBegin;
192   if (level > fas->levels-1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
193   if (fas->level !=  fas->levels - 1) SETERRQ2(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level);
194 
195   *lsnes = snes;
196   for (i = fas->level; i > level; i--) {
197     *lsnes = fas->next;
198     fas    = (SNES_FAS*)(*lsnes)->data;
199   }
200   if (fas->level != level) SETERRQ(PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
201   PetscFunctionReturn(0);
202 }
203 
204 #undef __FUNCT__
205 #define __FUNCT__ "SNESFASSetNumberSmoothUp"
206 /*@
207    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
208    use on all levels.
209 
210    Logically Collective on SNES
211 
212    Input Parameters:
213 +  snes - the multigrid context
214 -  n    - the number of smoothing steps
215 
216    Options Database Key:
217 .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
218 
219    Level: advanced
220 
221 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
222 
223 .seealso: SNESFASSetNumberSmoothDown()
224 @*/
225 PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
226 {
227   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
228   PetscErrorCode ierr;
229 
230   PetscFunctionBegin;
231   fas->max_up_it = n;
232   if (!fas->smoothu && fas->level != 0) {
233     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
234   }
235   if (fas->smoothu) {
236     ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr);
237   }
238   if (fas->next) {
239     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
240   }
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "SNESFASSetNumberSmoothDown"
246 /*@
247    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
248    use on all levels.
249 
250    Logically Collective on SNES
251 
252    Input Parameters:
253 +  snes - the multigrid context
254 -  n    - the number of smoothing steps
255 
256    Options Database Key:
257 .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
258 
259    Level: advanced
260 
261 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
262 
263 .seealso: SNESFASSetNumberSmoothUp()
264 @*/
265 PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
266 {
267   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
268   PetscErrorCode ierr = 0;
269 
270   PetscFunctionBegin;
271   if (!fas->smoothd) {
272     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
273   }
274   ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr);
275 
276   fas->max_down_it = n;
277   if (fas->next) {
278     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
279   }
280   PetscFunctionReturn(0);
281 }
282 
283 
284 #undef __FUNCT__
285 #define __FUNCT__ "SNESFASSetContinuation"
286 /*@
287    SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep
288 
289    Logically Collective on SNES
290 
291    Input Parameters:
292 +  snes - the multigrid context
293 -  n    - the number of smoothing steps
294 
295    Options Database Key:
296 .  -snes_fas_continuation - sets continuation to true
297 
298    Level: advanced
299 
300    Notes: This sets the prefix on the upsweep smoothers to -fas_continuation
301 
302 .keywords: FAS, MG, smoother, continuation
303 
304 .seealso: SNESFAS
305 @*/
306 PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation)
307 {
308   const char     *optionsprefix;
309   char           tprefix[128];
310   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
311   PetscErrorCode ierr = 0;
312 
313   PetscFunctionBegin;
314   ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr);
315   if (!fas->smoothu) {
316     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
317   }
318   sprintf(tprefix,"fas_levels_continuation_");
319   ierr = SNESSetOptionsPrefix(fas->smoothu, optionsprefix);CHKERRQ(ierr);
320   ierr = SNESAppendOptionsPrefix(fas->smoothu, tprefix);CHKERRQ(ierr);
321   ierr = SNESSetType(fas->smoothu,SNESNEWTONLS);CHKERRQ(ierr);
322   ierr = SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100);CHKERRQ(ierr);
323   fas->continuation = continuation;
324   if (fas->next) {
325     ierr = SNESFASSetContinuation(fas->next,continuation);CHKERRQ(ierr);
326   }
327   PetscFunctionReturn(0);
328 }
329 
330 
331 #undef __FUNCT__
332 #define __FUNCT__ "SNESFASSetCycles"
333 /*@
334    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
335    complicated cycling.
336 
337    Logically Collective on SNES
338 
339    Input Parameters:
340 +  snes   - the multigrid context
341 -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
342 
343    Options Database Key:
344 .  -snes_fas_cycles 1 or 2
345 
346    Level: advanced
347 
348 .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
349 
350 .seealso: SNESFASSetCyclesOnLevel()
351 @*/
352 PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
353 {
354   SNES_FAS       *fas = (SNES_FAS*)snes->data;
355   PetscErrorCode ierr;
356   PetscBool      isFine;
357 
358   PetscFunctionBegin;
359   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
360 
361   fas->n_cycles = cycles;
362   if (!isFine) {
363     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
364   }
365   if (fas->next) {
366     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
367   }
368   PetscFunctionReturn(0);
369 }
370 
371 
372 #undef __FUNCT__
373 #define __FUNCT__ "SNESFASSetMonitor"
374 /*@
375    SNESFASSetMonitor - Sets the method-specific cycle monitoring
376 
377    Logically Collective on SNES
378 
379    Input Parameters:
380 +  snes   - the FAS context
381 -  flg    - monitor or not
382 
383    Level: advanced
384 
385 .keywords: FAS, monitor
386 
387 .seealso: SNESFASSetCyclesOnLevel()
388 @*/
389 PetscErrorCode SNESFASSetMonitor(SNES snes, PetscBool flg)
390 {
391   SNES_FAS       *fas = (SNES_FAS*)snes->data;
392   PetscErrorCode ierr;
393   PetscBool      isFine;
394   PetscInt       i, levels = fas->levels;
395   SNES           levelsnes;
396 
397   PetscFunctionBegin;
398   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
399   if (isFine) {
400     for (i = 0; i < levels; i++) {
401       ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
402       fas  = (SNES_FAS*)levelsnes->data;
403       if (flg) {
404         fas->monitor = PETSC_VIEWER_STDOUT_(PetscObjectComm((PetscObject)levelsnes));CHKERRQ(ierr);
405         /* set the monitors for the upsmoother and downsmoother */
406         ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr);
407         ierr = SNESMonitorSet(levelsnes,SNESMonitorDefault,NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
408       } else if (i != fas->levels - 1) {
409         /* unset the monitors on the coarse levels */
410         ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr);
411       }
412     }
413   }
414   PetscFunctionReturn(0);
415 }
416 
417 #undef __FUNCT__
418 #define __FUNCT__ "SNESFASSetLog"
419 /*@
420    SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels
421 
422    Logically Collective on SNES
423 
424    Input Parameters:
425 +  snes   - the FAS context
426 -  flg    - monitor or not
427 
428    Level: advanced
429 
430 .keywords: FAS, logging
431 
432 .seealso: SNESFASSetMonitor()
433 @*/
434 PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
435 {
436   SNES_FAS       *fas = (SNES_FAS*)snes->data;
437   PetscErrorCode ierr;
438   PetscBool      isFine;
439   PetscInt       i, levels = fas->levels;
440   SNES           levelsnes;
441   char           eventname[128];
442 
443   PetscFunctionBegin;
444   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
445   if (isFine) {
446     for (i = 0; i < levels; i++) {
447       ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
448       fas  = (SNES_FAS*)levelsnes->data;
449       if (flg) {
450         sprintf(eventname,"FASSetup %d",(int)i);
451         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr);
452         sprintf(eventname,"FASSmooth %d",(int)i);
453         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr);
454         sprintf(eventname,"FASResid %d",(int)i);
455         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr);
456         if (i) {
457           sprintf(eventname,"FASInterp %d",(int)i);
458           ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr);
459         }
460       } else {
461         fas->eventsmoothsetup    = 0;
462         fas->eventsmoothsolve    = 0;
463         fas->eventresidual       = 0;
464         fas->eventinterprestrict = 0;
465       }
466     }
467   }
468   PetscFunctionReturn(0);
469 }
470 
471 #undef __FUNCT__
472 #define __FUNCT__ "SNESFASCycleCreateSmoother_Private"
473 /*
474 Creates the default smoother type.
475 
476 This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
477 
478  */
479 PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
480 {
481   SNES_FAS       *fas;
482   const char     *optionsprefix;
483   char           tprefix[128];
484   PetscErrorCode ierr;
485   SNES           nsmooth;
486 
487   PetscFunctionBegin;
488   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
489   fas  = (SNES_FAS*)snes->data;
490   ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr);
491   /* create the default smoother */
492   ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);CHKERRQ(ierr);
493   if (fas->level == 0) {
494     sprintf(tprefix,"fas_coarse_");
495     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
496     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
497     ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr);
498     ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr);
499   } else {
500     sprintf(tprefix,"fas_levels_%d_",(int)fas->level);
501     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
502     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
503     ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr);
504     ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr);
505   }
506   ierr    = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
507   ierr    = PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);CHKERRQ(ierr);
508   ierr    = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr);
509   ierr    = PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level);CHKERRQ(ierr);
510   *smooth = nsmooth;
511   PetscFunctionReturn(0);
512 }
513 
514 /* ------------- Functions called on a particular level ----------------- */
515 
516 #undef __FUNCT__
517 #define __FUNCT__ "SNESFASCycleSetCycles"
518 /*@
519    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
520 
521    Logically Collective on SNES
522 
523    Input Parameters:
524 +  snes   - the multigrid context
525 .  level  - the level to set the number of cycles on
526 -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
527 
528    Level: advanced
529 
530 .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid
531 
532 .seealso: SNESFASSetCycles()
533 @*/
534 PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
535 {
536   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
537   PetscErrorCode ierr;
538 
539   PetscFunctionBegin;
540   fas->n_cycles = cycles;
541   ierr          = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
542   PetscFunctionReturn(0);
543 }
544 
545 
546 #undef __FUNCT__
547 #define __FUNCT__ "SNESFASCycleGetSmoother"
548 /*@
549    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
550 
551    Logically Collective on SNES
552 
553    Input Parameters:
554 .  snes   - the multigrid context
555 
556    Output Parameters:
557 .  smooth - the smoother
558 
559    Level: advanced
560 
561 .keywords: SNES, FAS, get, smoother, multigrid
562 
563 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
564 @*/
565 PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
566 {
567   SNES_FAS *fas;
568 
569   PetscFunctionBegin;
570   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
571   fas     = (SNES_FAS*)snes->data;
572   *smooth = fas->smoothd;
573   PetscFunctionReturn(0);
574 }
575 #undef __FUNCT__
576 #define __FUNCT__ "SNESFASCycleGetSmootherUp"
577 /*@
578    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
579 
580    Logically Collective on SNES
581 
582    Input Parameters:
583 .  snes   - the multigrid context
584 
585    Output Parameters:
586 .  smoothu - the smoother
587 
588    Notes:
589    Returns the downsmoother if no up smoother is available.  This enables transparent
590    default behavior in the process of the solve.
591 
592    Level: advanced
593 
594 .keywords: SNES, FAS, get, smoother, multigrid
595 
596 .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
597 @*/
598 PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
599 {
600   SNES_FAS *fas;
601 
602   PetscFunctionBegin;
603   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
604   fas = (SNES_FAS*)snes->data;
605   if (!fas->smoothu) *smoothu = fas->smoothd;
606   else *smoothu = fas->smoothu;
607   PetscFunctionReturn(0);
608 }
609 
610 #undef __FUNCT__
611 #define __FUNCT__ "SNESFASCycleGetSmootherDown"
612 /*@
613    SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
614 
615    Logically Collective on SNES
616 
617    Input Parameters:
618 .  snes   - the multigrid context
619 
620    Output Parameters:
621 .  smoothd - the smoother
622 
623    Level: advanced
624 
625 .keywords: SNES, FAS, get, smoother, multigrid
626 
627 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
628 @*/
629 PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
630 {
631   SNES_FAS *fas;
632 
633   PetscFunctionBegin;
634   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
635   fas      = (SNES_FAS*)snes->data;
636   *smoothd = fas->smoothd;
637   PetscFunctionReturn(0);
638 }
639 
640 
641 #undef __FUNCT__
642 #define __FUNCT__ "SNESFASCycleGetCorrection"
643 /*@
644    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level
645 
646    Logically Collective on SNES
647 
648    Input Parameters:
649 .  snes   - the multigrid context
650 
651    Output Parameters:
652 .  correction - the coarse correction on this level
653 
654    Notes:
655    Returns NULL on the coarsest level.
656 
657    Level: advanced
658 
659 .keywords: SNES, FAS, get, smoother, multigrid
660 
661 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
662 @*/
663 PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
664 {
665   SNES_FAS *fas;
666 
667   PetscFunctionBegin;
668   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
669   fas         = (SNES_FAS*)snes->data;
670   *correction = fas->next;
671   PetscFunctionReturn(0);
672 }
673 
674 #undef __FUNCT__
675 #define __FUNCT__ "SNESFASCycleGetInterpolation"
676 /*@
677    SNESFASCycleGetInterpolation - Gets the interpolation on this level
678 
679    Logically Collective on SNES
680 
681    Input Parameters:
682 .  snes   - the multigrid context
683 
684    Output Parameters:
685 .  mat    - the interpolation operator on this level
686 
687    Level: developer
688 
689 .keywords: SNES, FAS, get, smoother, multigrid
690 
691 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
692 @*/
693 PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
694 {
695   SNES_FAS *fas;
696 
697   PetscFunctionBegin;
698   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
699   fas  = (SNES_FAS*)snes->data;
700   *mat = fas->interpolate;
701   PetscFunctionReturn(0);
702 }
703 
704 
705 #undef __FUNCT__
706 #define __FUNCT__ "SNESFASCycleGetRestriction"
707 /*@
708    SNESFASCycleGetRestriction - Gets the restriction on this level
709 
710    Logically Collective on SNES
711 
712    Input Parameters:
713 .  snes   - the multigrid context
714 
715    Output Parameters:
716 .  mat    - the restriction operator on this level
717 
718    Level: developer
719 
720 .keywords: SNES, FAS, get, smoother, multigrid
721 
722 .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
723 @*/
724 PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
725 {
726   SNES_FAS *fas;
727 
728   PetscFunctionBegin;
729   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
730   fas  = (SNES_FAS*)snes->data;
731   *mat = fas->restrct;
732   PetscFunctionReturn(0);
733 }
734 
735 
736 #undef __FUNCT__
737 #define __FUNCT__ "SNESFASCycleGetInjection"
738 /*@
739    SNESFASCycleGetInjection - Gets the injection on this level
740 
741    Logically Collective on SNES
742 
743    Input Parameters:
744 .  snes   - the multigrid context
745 
746    Output Parameters:
747 .  mat    - the restriction operator on this level
748 
749    Level: developer
750 
751 .keywords: SNES, FAS, get, smoother, multigrid
752 
753 .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
754 @*/
755 PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
756 {
757   SNES_FAS *fas;
758 
759   PetscFunctionBegin;
760   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
761   fas  = (SNES_FAS*)snes->data;
762   *mat = fas->inject;
763   PetscFunctionReturn(0);
764 }
765 
766 #undef __FUNCT__
767 #define __FUNCT__ "SNESFASCycleGetRScale"
768 /*@
769    SNESFASCycleGetRScale - Gets the injection on this level
770 
771    Logically Collective on SNES
772 
773    Input Parameters:
774 .  snes   - the multigrid context
775 
776    Output Parameters:
777 .  mat    - the restriction operator on this level
778 
779    Level: developer
780 
781 .keywords: SNES, FAS, get, smoother, multigrid
782 
783 .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
784 @*/
785 PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
786 {
787   SNES_FAS *fas;
788 
789   PetscFunctionBegin;
790   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
791   fas  = (SNES_FAS*)snes->data;
792   *vec = fas->rscale;
793   PetscFunctionReturn(0);
794 }
795 
796 #undef __FUNCT__
797 #define __FUNCT__ "SNESFASCycleIsFine"
798 /*@
799    SNESFASCycleIsFine - Determines if a given cycle is the fine level.
800 
801    Logically Collective on SNES
802 
803    Input Parameters:
804 .  snes   - the FAS context
805 
806    Output Parameters:
807 .  flg - indicates if this is the fine level or not
808 
809    Level: advanced
810 
811 .keywords: SNES, FAS
812 
813 .seealso: SNESFASSetLevels()
814 @*/
815 PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
816 {
817   SNES_FAS *fas;
818 
819   PetscFunctionBegin;
820   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
821   fas = (SNES_FAS*)snes->data;
822   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
823   else *flg = PETSC_FALSE;
824   PetscFunctionReturn(0);
825 }
826 
827 /* ---------- functions called on the finest level that return level-specific information ---------- */
828 
829 #undef __FUNCT__
830 #define __FUNCT__ "SNESFASSetInterpolation"
831 /*@
832    SNESFASSetInterpolation - Sets the function to be used to calculate the
833    interpolation from l-1 to the lth level
834 
835    Input Parameters:
836 +  snes      - the multigrid context
837 .  mat       - the interpolation operator
838 -  level     - the level (0 is coarsest) to supply [do not supply 0]
839 
840    Level: advanced
841 
842    Notes:
843           Usually this is the same matrix used also to set the restriction
844     for the same level.
845 
846           One can pass in the interpolation matrix or its transpose; PETSc figures
847     out from the matrix size which one it is.
848 
849 .keywords:  FAS, multigrid, set, interpolate, level
850 
851 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
852 @*/
853 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
854 {
855   SNES_FAS       *fas;
856   PetscErrorCode ierr;
857   SNES           levelsnes;
858 
859   PetscFunctionBegin;
860   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
861   fas  = (SNES_FAS*)levelsnes->data;
862   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
863   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
864 
865   fas->interpolate = mat;
866   PetscFunctionReturn(0);
867 }
868 
869 #undef __FUNCT__
870 #define __FUNCT__ "SNESFASGetInterpolation"
871 /*@
872    SNESFASGetInterpolation - Gets the matrix used to calculate the
873    interpolation from l-1 to the lth level
874 
875    Input Parameters:
876 +  snes      - the multigrid context
877 -  level     - the level (0 is coarsest) to supply [do not supply 0]
878 
879    Output Parameters:
880 .  mat       - the interpolation operator
881 
882    Level: advanced
883 
884 .keywords:  FAS, multigrid, get, interpolate, level
885 
886 .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
887 @*/
888 PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
889 {
890   SNES_FAS       *fas;
891   PetscErrorCode ierr;
892   SNES           levelsnes;
893 
894   PetscFunctionBegin;
895   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
896   fas  = (SNES_FAS*)levelsnes->data;
897   *mat = fas->interpolate;
898   PetscFunctionReturn(0);
899 }
900 
901 #undef __FUNCT__
902 #define __FUNCT__ "SNESFASSetRestriction"
903 /*@
904    SNESFASSetRestriction - Sets the function to be used to restrict the defect
905    from level l to l-1.
906 
907    Input Parameters:
908 +  snes  - the multigrid context
909 .  mat   - the restriction matrix
910 -  level - the level (0 is coarsest) to supply [Do not supply 0]
911 
912    Level: advanced
913 
914    Notes:
915           Usually this is the same matrix used also to set the interpolation
916     for the same level.
917 
918           One can pass in the interpolation matrix or its transpose; PETSc figures
919     out from the matrix size which one it is.
920 
921          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
922     is used.
923 
924 .keywords: FAS, MG, set, multigrid, restriction, level
925 
926 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
927 @*/
928 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
929 {
930   SNES_FAS       *fas;
931   PetscErrorCode ierr;
932   SNES           levelsnes;
933 
934   PetscFunctionBegin;
935   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
936   fas  = (SNES_FAS*)levelsnes->data;
937   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
938   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
939 
940   fas->restrct = mat;
941   PetscFunctionReturn(0);
942 }
943 
944 #undef __FUNCT__
945 #define __FUNCT__ "SNESFASGetRestriction"
946 /*@
947    SNESFASGetRestriction - Gets the matrix used to calculate the
948    restriction from l to the l-1th level
949 
950    Input Parameters:
951 +  snes      - the multigrid context
952 -  level     - the level (0 is coarsest) to supply [do not supply 0]
953 
954    Output Parameters:
955 .  mat       - the interpolation operator
956 
957    Level: advanced
958 
959 .keywords:  FAS, multigrid, get, restrict, level
960 
961 .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
962 @*/
963 PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
964 {
965   SNES_FAS       *fas;
966   PetscErrorCode ierr;
967   SNES           levelsnes;
968 
969   PetscFunctionBegin;
970   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
971   fas  = (SNES_FAS*)levelsnes->data;
972   *mat = fas->restrct;
973   PetscFunctionReturn(0);
974 }
975 
976 
977 #undef __FUNCT__
978 #define __FUNCT__ "SNESFASSetInjection"
979 /*@
980    SNESFASSetInjection - Sets the function to be used to inject the solution
981    from level l to l-1.
982 
983    Input Parameters:
984  +  snes  - the multigrid context
985 .  mat   - the restriction matrix
986 -  level - the level (0 is coarsest) to supply [Do not supply 0]
987 
988    Level: advanced
989 
990    Notes:
991          If you do not set this, the restriction and rscale is used to
992    project the solution instead.
993 
994 .keywords: FAS, MG, set, multigrid, restriction, level
995 
996 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
997 @*/
998 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
999 {
1000   SNES_FAS       *fas;
1001   PetscErrorCode ierr;
1002   SNES           levelsnes;
1003 
1004   PetscFunctionBegin;
1005   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1006   fas  = (SNES_FAS*)levelsnes->data;
1007   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1008   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
1009 
1010   fas->inject = mat;
1011   PetscFunctionReturn(0);
1012 }
1013 
1014 
1015 #undef __FUNCT__
1016 #define __FUNCT__ "SNESFASGetInjection"
1017 /*@
1018    SNESFASGetInjection - Gets the matrix used to calculate the
1019    injection from l-1 to the lth level
1020 
1021    Input Parameters:
1022 +  snes      - the multigrid context
1023 -  level     - the level (0 is coarsest) to supply [do not supply 0]
1024 
1025    Output Parameters:
1026 .  mat       - the injection operator
1027 
1028    Level: advanced
1029 
1030 .keywords:  FAS, multigrid, get, restrict, level
1031 
1032 .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
1033 @*/
1034 PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
1035 {
1036   SNES_FAS       *fas;
1037   PetscErrorCode ierr;
1038   SNES           levelsnes;
1039 
1040   PetscFunctionBegin;
1041   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1042   fas  = (SNES_FAS*)levelsnes->data;
1043   *mat = fas->inject;
1044   PetscFunctionReturn(0);
1045 }
1046 
1047 #undef __FUNCT__
1048 #define __FUNCT__ "SNESFASSetRScale"
1049 /*@
1050    SNESFASSetRScale - Sets the scaling factor of the restriction
1051    operator from level l to l-1.
1052 
1053    Input Parameters:
1054 +  snes   - the multigrid context
1055 .  rscale - the restriction scaling
1056 -  level  - the level (0 is coarsest) to supply [Do not supply 0]
1057 
1058    Level: advanced
1059 
1060    Notes:
1061          This is only used in the case that the injection is not set.
1062 
1063 .keywords: FAS, MG, set, multigrid, restriction, level
1064 
1065 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1066 @*/
1067 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
1068 {
1069   SNES_FAS       *fas;
1070   PetscErrorCode ierr;
1071   SNES           levelsnes;
1072 
1073   PetscFunctionBegin;
1074   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1075   fas  = (SNES_FAS*)levelsnes->data;
1076   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
1077   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
1078 
1079   fas->rscale = rscale;
1080   PetscFunctionReturn(0);
1081 }
1082 
1083 #undef __FUNCT__
1084 #define __FUNCT__ "SNESFASGetSmoother"
1085 /*@
1086    SNESFASGetSmoother - Gets the default smoother on a level.
1087 
1088    Input Parameters:
1089 +  snes   - the multigrid context
1090 -  level  - the level (0 is coarsest) to supply
1091 
1092    Output Parameters:
1093    smooth  - the smoother
1094 
1095    Level: advanced
1096 
1097 .keywords: FAS, MG, get, multigrid, smoother, level
1098 
1099 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1100 @*/
1101 PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
1102 {
1103   SNES_FAS       *fas;
1104   PetscErrorCode ierr;
1105   SNES           levelsnes;
1106 
1107   PetscFunctionBegin;
1108   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1109   fas  = (SNES_FAS*)levelsnes->data;
1110   if (!fas->smoothd) {
1111     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr);
1112   }
1113   *smooth = fas->smoothd;
1114   PetscFunctionReturn(0);
1115 }
1116 
1117 #undef __FUNCT__
1118 #define __FUNCT__ "SNESFASGetSmootherDown"
1119 /*@
1120    SNESFASGetSmootherDown - Gets the downsmoother on a level.
1121 
1122    Input Parameters:
1123 +  snes   - the multigrid context
1124 -  level  - the level (0 is coarsest) to supply
1125 
1126    Output Parameters:
1127    smooth  - the smoother
1128 
1129    Level: advanced
1130 
1131 .keywords: FAS, MG, get, multigrid, smoother, level
1132 
1133 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1134 @*/
1135 PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
1136 {
1137   SNES_FAS       *fas;
1138   PetscErrorCode ierr;
1139   SNES           levelsnes;
1140 
1141   PetscFunctionBegin;
1142   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1143   fas  = (SNES_FAS*)levelsnes->data;
1144   /* if the user chooses to differentiate smoothers, create them both at this point */
1145   if (!fas->smoothd) {
1146     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr);
1147   }
1148   if (!fas->smoothu) {
1149     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr);
1150   }
1151   *smooth = fas->smoothd;
1152   PetscFunctionReturn(0);
1153 }
1154 
1155 #undef __FUNCT__
1156 #define __FUNCT__ "SNESFASGetSmootherUp"
1157 /*@
1158    SNESFASGetSmootherUp - Gets the upsmoother on a level.
1159 
1160    Input Parameters:
1161 +  snes   - the multigrid context
1162 -  level  - the level (0 is coarsest)
1163 
1164    Output Parameters:
1165    smooth  - the smoother
1166 
1167    Level: advanced
1168 
1169 .keywords: FAS, MG, get, multigrid, smoother, level
1170 
1171 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1172 @*/
1173 PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
1174 {
1175   SNES_FAS       *fas;
1176   PetscErrorCode ierr;
1177   SNES           levelsnes;
1178 
1179   PetscFunctionBegin;
1180   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1181   fas  = (SNES_FAS*)levelsnes->data;
1182   /* if the user chooses to differentiate smoothers, create them both at this point */
1183   if (!fas->smoothd) {
1184     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr);
1185   }
1186   if (!fas->smoothu) {
1187     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr);
1188   }
1189   *smooth = fas->smoothu;
1190   PetscFunctionReturn(0);
1191 }
1192 
1193 #undef __FUNCT__
1194 #define __FUNCT__ "SNESFASGetCoarseSolve"
1195 /*@
1196   SNESFASGetCoarseSolve - Gets the coarsest solver.
1197 
1198   Input Parameters:
1199 . snes - the multigrid context
1200 
1201   Output Parameters:
1202 . coarse - the coarse-level solver
1203 
1204   Level: advanced
1205 
1206 .keywords: FAS, MG, get, multigrid, solver, coarse
1207 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1208 @*/
1209 PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
1210 {
1211   SNES_FAS       *fas;
1212   PetscErrorCode ierr;
1213   SNES           levelsnes;
1214 
1215   PetscFunctionBegin;
1216   ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr);
1217   fas  = (SNES_FAS*)levelsnes->data;
1218   /* if the user chooses to differentiate smoothers, create them both at this point */
1219   if (!fas->smoothd) {
1220     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr);
1221   }
1222   *coarse = fas->smoothd;
1223   PetscFunctionReturn(0);
1224 }
1225 
1226 #undef __FUNCT__
1227 #define __FUNCT__ "SNESFASFullSetDownSweep"
1228 /*@
1229    SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS
1230 
1231    Logically Collective on SNES
1232 
1233    Input Parameters:
1234 +  snes - the multigrid context
1235 -  swp - whether to downsweep or not
1236 
1237    Options Database Key:
1238 .  -snes_fas_full_downsweep - Sets number of pre-smoothing steps
1239 
1240    Level: advanced
1241 
1242 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
1243 
1244 .seealso: SNESFASSetNumberSmoothUp()
1245 @*/
1246 PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp)
1247 {
1248   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
1249   PetscErrorCode ierr = 0;
1250 
1251   PetscFunctionBegin;
1252   fas->full_downsweep = swp;
1253   if (fas->next) {
1254     ierr = SNESFASFullSetDownSweep(fas->next,swp);CHKERRQ(ierr);
1255   }
1256   PetscFunctionReturn(0);
1257 }
1258