xref: /petsc/src/snes/impls/fas/fasfunc.c (revision 70baa948e42e6bb1660071283af10f413c9269b6)
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 .  vf     - viewer and format structure (may be NULL if flg is FALSE)
382 -  flg    - monitor or not
383 
384    Level: advanced
385 
386 .keywords: FAS, monitor
387 
388 .seealso: SNESFASSetCyclesOnLevel()
389 @*/
390 PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg)
391 {
392   SNES_FAS       *fas = (SNES_FAS*)snes->data;
393   PetscErrorCode ierr;
394   PetscBool      isFine;
395   PetscInt       i, levels = fas->levels;
396   SNES           levelsnes;
397 
398   PetscFunctionBegin;
399   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
400   if (isFine) {
401     for (i = 0; i < levels; i++) {
402       ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
403       fas  = (SNES_FAS*)levelsnes->data;
404       if (flg) {
405         /* set the monitors for the upsmoother and downsmoother */
406         ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr);
407         /* Only register destroy on finest level */
408         ierr = SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL));CHKERRQ(ierr);
409       } else if (i != fas->levels - 1) {
410         /* unset the monitors on the coarse levels */
411         ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr);
412       }
413     }
414   }
415   PetscFunctionReturn(0);
416 }
417 
418 #undef __FUNCT__
419 #define __FUNCT__ "SNESFASSetLog"
420 /*@
421    SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels
422 
423    Logically Collective on SNES
424 
425    Input Parameters:
426 +  snes   - the FAS context
427 -  flg    - monitor or not
428 
429    Level: advanced
430 
431 .keywords: FAS, logging
432 
433 .seealso: SNESFASSetMonitor()
434 @*/
435 PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
436 {
437   SNES_FAS       *fas = (SNES_FAS*)snes->data;
438   PetscErrorCode ierr;
439   PetscBool      isFine;
440   PetscInt       i, levels = fas->levels;
441   SNES           levelsnes;
442   char           eventname[128];
443 
444   PetscFunctionBegin;
445   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
446   if (isFine) {
447     for (i = 0; i < levels; i++) {
448       ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
449       fas  = (SNES_FAS*)levelsnes->data;
450       if (flg) {
451         sprintf(eventname,"FASSetup %d",(int)i);
452         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr);
453         sprintf(eventname,"FASSmooth %d",(int)i);
454         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr);
455         sprintf(eventname,"FASResid %d",(int)i);
456         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr);
457         if (i) {
458           sprintf(eventname,"FASInterp %d",(int)i);
459           ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr);
460         }
461       } else {
462         fas->eventsmoothsetup    = 0;
463         fas->eventsmoothsolve    = 0;
464         fas->eventresidual       = 0;
465         fas->eventinterprestrict = 0;
466       }
467     }
468   }
469   PetscFunctionReturn(0);
470 }
471 
472 #undef __FUNCT__
473 #define __FUNCT__ "SNESFASCycleCreateSmoother_Private"
474 /*
475 Creates the default smoother type.
476 
477 This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
478 
479  */
480 PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
481 {
482   SNES_FAS       *fas;
483   const char     *optionsprefix;
484   char           tprefix[128];
485   PetscErrorCode ierr;
486   SNES           nsmooth;
487 
488   PetscFunctionBegin;
489   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
490   fas  = (SNES_FAS*)snes->data;
491   ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr);
492   /* create the default smoother */
493   ierr = SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth);CHKERRQ(ierr);
494   if (fas->level == 0) {
495     sprintf(tprefix,"fas_coarse_");
496     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
497     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
498     ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr);
499     ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr);
500   } else {
501     sprintf(tprefix,"fas_levels_%d_",(int)fas->level);
502     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
503     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
504     ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr);
505     ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr);
506   }
507   ierr    = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
508   ierr    = PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth);CHKERRQ(ierr);
509   ierr    = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr);
510   ierr    = PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level);CHKERRQ(ierr);
511   *smooth = nsmooth;
512   PetscFunctionReturn(0);
513 }
514 
515 /* ------------- Functions called on a particular level ----------------- */
516 
517 #undef __FUNCT__
518 #define __FUNCT__ "SNESFASCycleSetCycles"
519 /*@
520    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
521 
522    Logically Collective on SNES
523 
524    Input Parameters:
525 +  snes   - the multigrid context
526 .  level  - the level to set the number of cycles on
527 -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
528 
529    Level: advanced
530 
531 .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid
532 
533 .seealso: SNESFASSetCycles()
534 @*/
535 PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
536 {
537   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
538   PetscErrorCode ierr;
539 
540   PetscFunctionBegin;
541   fas->n_cycles = cycles;
542   ierr          = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
543   PetscFunctionReturn(0);
544 }
545 
546 
547 #undef __FUNCT__
548 #define __FUNCT__ "SNESFASCycleGetSmoother"
549 /*@
550    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
551 
552    Logically Collective on SNES
553 
554    Input Parameters:
555 .  snes   - the multigrid context
556 
557    Output Parameters:
558 .  smooth - the smoother
559 
560    Level: advanced
561 
562 .keywords: SNES, FAS, get, smoother, multigrid
563 
564 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
565 @*/
566 PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
567 {
568   SNES_FAS *fas;
569 
570   PetscFunctionBegin;
571   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
572   fas     = (SNES_FAS*)snes->data;
573   *smooth = fas->smoothd;
574   PetscFunctionReturn(0);
575 }
576 #undef __FUNCT__
577 #define __FUNCT__ "SNESFASCycleGetSmootherUp"
578 /*@
579    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
580 
581    Logically Collective on SNES
582 
583    Input Parameters:
584 .  snes   - the multigrid context
585 
586    Output Parameters:
587 .  smoothu - the smoother
588 
589    Notes:
590    Returns the downsmoother if no up smoother is available.  This enables transparent
591    default behavior in the process of the solve.
592 
593    Level: advanced
594 
595 .keywords: SNES, FAS, get, smoother, multigrid
596 
597 .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
598 @*/
599 PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
600 {
601   SNES_FAS *fas;
602 
603   PetscFunctionBegin;
604   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
605   fas = (SNES_FAS*)snes->data;
606   if (!fas->smoothu) *smoothu = fas->smoothd;
607   else *smoothu = fas->smoothu;
608   PetscFunctionReturn(0);
609 }
610 
611 #undef __FUNCT__
612 #define __FUNCT__ "SNESFASCycleGetSmootherDown"
613 /*@
614    SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
615 
616    Logically Collective on SNES
617 
618    Input Parameters:
619 .  snes   - the multigrid context
620 
621    Output Parameters:
622 .  smoothd - the smoother
623 
624    Level: advanced
625 
626 .keywords: SNES, FAS, get, smoother, multigrid
627 
628 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
629 @*/
630 PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
631 {
632   SNES_FAS *fas;
633 
634   PetscFunctionBegin;
635   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
636   fas      = (SNES_FAS*)snes->data;
637   *smoothd = fas->smoothd;
638   PetscFunctionReturn(0);
639 }
640 
641 
642 #undef __FUNCT__
643 #define __FUNCT__ "SNESFASCycleGetCorrection"
644 /*@
645    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level
646 
647    Logically Collective on SNES
648 
649    Input Parameters:
650 .  snes   - the multigrid context
651 
652    Output Parameters:
653 .  correction - the coarse correction on this level
654 
655    Notes:
656    Returns NULL on the coarsest level.
657 
658    Level: advanced
659 
660 .keywords: SNES, FAS, get, smoother, multigrid
661 
662 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
663 @*/
664 PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
665 {
666   SNES_FAS *fas;
667 
668   PetscFunctionBegin;
669   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
670   fas         = (SNES_FAS*)snes->data;
671   *correction = fas->next;
672   PetscFunctionReturn(0);
673 }
674 
675 #undef __FUNCT__
676 #define __FUNCT__ "SNESFASCycleGetInterpolation"
677 /*@
678    SNESFASCycleGetInterpolation - Gets the interpolation on this level
679 
680    Logically Collective on SNES
681 
682    Input Parameters:
683 .  snes   - the multigrid context
684 
685    Output Parameters:
686 .  mat    - the interpolation operator on this level
687 
688    Level: developer
689 
690 .keywords: SNES, FAS, get, smoother, multigrid
691 
692 .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
693 @*/
694 PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
695 {
696   SNES_FAS *fas;
697 
698   PetscFunctionBegin;
699   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
700   fas  = (SNES_FAS*)snes->data;
701   *mat = fas->interpolate;
702   PetscFunctionReturn(0);
703 }
704 
705 
706 #undef __FUNCT__
707 #define __FUNCT__ "SNESFASCycleGetRestriction"
708 /*@
709    SNESFASCycleGetRestriction - Gets the restriction on this level
710 
711    Logically Collective on SNES
712 
713    Input Parameters:
714 .  snes   - the multigrid context
715 
716    Output Parameters:
717 .  mat    - the restriction operator on this level
718 
719    Level: developer
720 
721 .keywords: SNES, FAS, get, smoother, multigrid
722 
723 .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
724 @*/
725 PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
726 {
727   SNES_FAS *fas;
728 
729   PetscFunctionBegin;
730   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
731   fas  = (SNES_FAS*)snes->data;
732   *mat = fas->restrct;
733   PetscFunctionReturn(0);
734 }
735 
736 
737 #undef __FUNCT__
738 #define __FUNCT__ "SNESFASCycleGetInjection"
739 /*@
740    SNESFASCycleGetInjection - Gets the injection on this level
741 
742    Logically Collective on SNES
743 
744    Input Parameters:
745 .  snes   - the multigrid context
746 
747    Output Parameters:
748 .  mat    - the restriction operator on this level
749 
750    Level: developer
751 
752 .keywords: SNES, FAS, get, smoother, multigrid
753 
754 .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
755 @*/
756 PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
757 {
758   SNES_FAS *fas;
759 
760   PetscFunctionBegin;
761   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
762   fas  = (SNES_FAS*)snes->data;
763   *mat = fas->inject;
764   PetscFunctionReturn(0);
765 }
766 
767 #undef __FUNCT__
768 #define __FUNCT__ "SNESFASCycleGetRScale"
769 /*@
770    SNESFASCycleGetRScale - Gets the injection on this level
771 
772    Logically Collective on SNES
773 
774    Input Parameters:
775 .  snes   - the multigrid context
776 
777    Output Parameters:
778 .  mat    - the restriction operator on this level
779 
780    Level: developer
781 
782 .keywords: SNES, FAS, get, smoother, multigrid
783 
784 .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
785 @*/
786 PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
787 {
788   SNES_FAS *fas;
789 
790   PetscFunctionBegin;
791   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
792   fas  = (SNES_FAS*)snes->data;
793   *vec = fas->rscale;
794   PetscFunctionReturn(0);
795 }
796 
797 #undef __FUNCT__
798 #define __FUNCT__ "SNESFASCycleIsFine"
799 /*@
800    SNESFASCycleIsFine - Determines if a given cycle is the fine level.
801 
802    Logically Collective on SNES
803 
804    Input Parameters:
805 .  snes   - the FAS context
806 
807    Output Parameters:
808 .  flg - indicates if this is the fine level or not
809 
810    Level: advanced
811 
812 .keywords: SNES, FAS
813 
814 .seealso: SNESFASSetLevels()
815 @*/
816 PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
817 {
818   SNES_FAS *fas;
819 
820   PetscFunctionBegin;
821   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
822   fas = (SNES_FAS*)snes->data;
823   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
824   else *flg = PETSC_FALSE;
825   PetscFunctionReturn(0);
826 }
827 
828 /* ---------- functions called on the finest level that return level-specific information ---------- */
829 
830 #undef __FUNCT__
831 #define __FUNCT__ "SNESFASSetInterpolation"
832 /*@
833    SNESFASSetInterpolation - Sets the function to be used to calculate the
834    interpolation from l-1 to the lth level
835 
836    Input Parameters:
837 +  snes      - the multigrid context
838 .  mat       - the interpolation operator
839 -  level     - the level (0 is coarsest) to supply [do not supply 0]
840 
841    Level: advanced
842 
843    Notes:
844           Usually this is the same matrix used also to set the restriction
845     for the same level.
846 
847           One can pass in the interpolation matrix or its transpose; PETSc figures
848     out from the matrix size which one it is.
849 
850 .keywords:  FAS, multigrid, set, interpolate, level
851 
852 .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
853 @*/
854 PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
855 {
856   SNES_FAS       *fas;
857   PetscErrorCode ierr;
858   SNES           levelsnes;
859 
860   PetscFunctionBegin;
861   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
862   fas  = (SNES_FAS*)levelsnes->data;
863   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
864   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
865 
866   fas->interpolate = mat;
867   PetscFunctionReturn(0);
868 }
869 
870 #undef __FUNCT__
871 #define __FUNCT__ "SNESFASGetInterpolation"
872 /*@
873    SNESFASGetInterpolation - Gets the matrix used to calculate the
874    interpolation from l-1 to the lth level
875 
876    Input Parameters:
877 +  snes      - the multigrid context
878 -  level     - the level (0 is coarsest) to supply [do not supply 0]
879 
880    Output Parameters:
881 .  mat       - the interpolation operator
882 
883    Level: advanced
884 
885 .keywords:  FAS, multigrid, get, interpolate, level
886 
887 .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
888 @*/
889 PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
890 {
891   SNES_FAS       *fas;
892   PetscErrorCode ierr;
893   SNES           levelsnes;
894 
895   PetscFunctionBegin;
896   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
897   fas  = (SNES_FAS*)levelsnes->data;
898   *mat = fas->interpolate;
899   PetscFunctionReturn(0);
900 }
901 
902 #undef __FUNCT__
903 #define __FUNCT__ "SNESFASSetRestriction"
904 /*@
905    SNESFASSetRestriction - Sets the function to be used to restrict the defect
906    from level l to l-1.
907 
908    Input Parameters:
909 +  snes  - the multigrid context
910 .  mat   - the restriction matrix
911 -  level - the level (0 is coarsest) to supply [Do not supply 0]
912 
913    Level: advanced
914 
915    Notes:
916           Usually this is the same matrix used also to set the interpolation
917     for the same level.
918 
919           One can pass in the interpolation matrix or its transpose; PETSc figures
920     out from the matrix size which one it is.
921 
922          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
923     is used.
924 
925 .keywords: FAS, MG, set, multigrid, restriction, level
926 
927 .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
928 @*/
929 PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
930 {
931   SNES_FAS       *fas;
932   PetscErrorCode ierr;
933   SNES           levelsnes;
934 
935   PetscFunctionBegin;
936   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
937   fas  = (SNES_FAS*)levelsnes->data;
938   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
939   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
940 
941   fas->restrct = mat;
942   PetscFunctionReturn(0);
943 }
944 
945 #undef __FUNCT__
946 #define __FUNCT__ "SNESFASGetRestriction"
947 /*@
948    SNESFASGetRestriction - Gets the matrix used to calculate the
949    restriction from l to the l-1th level
950 
951    Input Parameters:
952 +  snes      - the multigrid context
953 -  level     - the level (0 is coarsest) to supply [do not supply 0]
954 
955    Output Parameters:
956 .  mat       - the interpolation operator
957 
958    Level: advanced
959 
960 .keywords:  FAS, multigrid, get, restrict, level
961 
962 .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
963 @*/
964 PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
965 {
966   SNES_FAS       *fas;
967   PetscErrorCode ierr;
968   SNES           levelsnes;
969 
970   PetscFunctionBegin;
971   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
972   fas  = (SNES_FAS*)levelsnes->data;
973   *mat = fas->restrct;
974   PetscFunctionReturn(0);
975 }
976 
977 
978 #undef __FUNCT__
979 #define __FUNCT__ "SNESFASSetInjection"
980 /*@
981    SNESFASSetInjection - Sets the function to be used to inject the solution
982    from level l to l-1.
983 
984    Input Parameters:
985  +  snes  - the multigrid context
986 .  mat   - the restriction matrix
987 -  level - the level (0 is coarsest) to supply [Do not supply 0]
988 
989    Level: advanced
990 
991    Notes:
992          If you do not set this, the restriction and rscale is used to
993    project the solution instead.
994 
995 .keywords: FAS, MG, set, multigrid, restriction, level
996 
997 .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
998 @*/
999 PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
1000 {
1001   SNES_FAS       *fas;
1002   PetscErrorCode ierr;
1003   SNES           levelsnes;
1004 
1005   PetscFunctionBegin;
1006   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1007   fas  = (SNES_FAS*)levelsnes->data;
1008   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
1009   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
1010 
1011   fas->inject = mat;
1012   PetscFunctionReturn(0);
1013 }
1014 
1015 
1016 #undef __FUNCT__
1017 #define __FUNCT__ "SNESFASGetInjection"
1018 /*@
1019    SNESFASGetInjection - Gets the matrix used to calculate the
1020    injection from l-1 to the lth level
1021 
1022    Input Parameters:
1023 +  snes      - the multigrid context
1024 -  level     - the level (0 is coarsest) to supply [do not supply 0]
1025 
1026    Output Parameters:
1027 .  mat       - the injection operator
1028 
1029    Level: advanced
1030 
1031 .keywords:  FAS, multigrid, get, restrict, level
1032 
1033 .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
1034 @*/
1035 PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
1036 {
1037   SNES_FAS       *fas;
1038   PetscErrorCode ierr;
1039   SNES           levelsnes;
1040 
1041   PetscFunctionBegin;
1042   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1043   fas  = (SNES_FAS*)levelsnes->data;
1044   *mat = fas->inject;
1045   PetscFunctionReturn(0);
1046 }
1047 
1048 #undef __FUNCT__
1049 #define __FUNCT__ "SNESFASSetRScale"
1050 /*@
1051    SNESFASSetRScale - Sets the scaling factor of the restriction
1052    operator from level l to l-1.
1053 
1054    Input Parameters:
1055 +  snes   - the multigrid context
1056 .  rscale - the restriction scaling
1057 -  level  - the level (0 is coarsest) to supply [Do not supply 0]
1058 
1059    Level: advanced
1060 
1061    Notes:
1062          This is only used in the case that the injection is not set.
1063 
1064 .keywords: FAS, MG, set, multigrid, restriction, level
1065 
1066 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1067 @*/
1068 PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
1069 {
1070   SNES_FAS       *fas;
1071   PetscErrorCode ierr;
1072   SNES           levelsnes;
1073 
1074   PetscFunctionBegin;
1075   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1076   fas  = (SNES_FAS*)levelsnes->data;
1077   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
1078   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
1079 
1080   fas->rscale = rscale;
1081   PetscFunctionReturn(0);
1082 }
1083 
1084 #undef __FUNCT__
1085 #define __FUNCT__ "SNESFASGetSmoother"
1086 /*@
1087    SNESFASGetSmoother - Gets the default smoother on a level.
1088 
1089    Input Parameters:
1090 +  snes   - the multigrid context
1091 -  level  - the level (0 is coarsest) to supply
1092 
1093    Output Parameters:
1094    smooth  - the smoother
1095 
1096    Level: advanced
1097 
1098 .keywords: FAS, MG, get, multigrid, smoother, level
1099 
1100 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1101 @*/
1102 PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
1103 {
1104   SNES_FAS       *fas;
1105   PetscErrorCode ierr;
1106   SNES           levelsnes;
1107 
1108   PetscFunctionBegin;
1109   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1110   fas  = (SNES_FAS*)levelsnes->data;
1111   if (!fas->smoothd) {
1112     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr);
1113   }
1114   *smooth = fas->smoothd;
1115   PetscFunctionReturn(0);
1116 }
1117 
1118 #undef __FUNCT__
1119 #define __FUNCT__ "SNESFASGetSmootherDown"
1120 /*@
1121    SNESFASGetSmootherDown - Gets the downsmoother on a level.
1122 
1123    Input Parameters:
1124 +  snes   - the multigrid context
1125 -  level  - the level (0 is coarsest) to supply
1126 
1127    Output Parameters:
1128    smooth  - the smoother
1129 
1130    Level: advanced
1131 
1132 .keywords: FAS, MG, get, multigrid, smoother, level
1133 
1134 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1135 @*/
1136 PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
1137 {
1138   SNES_FAS       *fas;
1139   PetscErrorCode ierr;
1140   SNES           levelsnes;
1141 
1142   PetscFunctionBegin;
1143   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1144   fas  = (SNES_FAS*)levelsnes->data;
1145   /* if the user chooses to differentiate smoothers, create them both at this point */
1146   if (!fas->smoothd) {
1147     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr);
1148   }
1149   if (!fas->smoothu) {
1150     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr);
1151   }
1152   *smooth = fas->smoothd;
1153   PetscFunctionReturn(0);
1154 }
1155 
1156 #undef __FUNCT__
1157 #define __FUNCT__ "SNESFASGetSmootherUp"
1158 /*@
1159    SNESFASGetSmootherUp - Gets the upsmoother on a level.
1160 
1161    Input Parameters:
1162 +  snes   - the multigrid context
1163 -  level  - the level (0 is coarsest)
1164 
1165    Output Parameters:
1166    smooth  - the smoother
1167 
1168    Level: advanced
1169 
1170 .keywords: FAS, MG, get, multigrid, smoother, level
1171 
1172 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1173 @*/
1174 PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
1175 {
1176   SNES_FAS       *fas;
1177   PetscErrorCode ierr;
1178   SNES           levelsnes;
1179 
1180   PetscFunctionBegin;
1181   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1182   fas  = (SNES_FAS*)levelsnes->data;
1183   /* if the user chooses to differentiate smoothers, create them both at this point */
1184   if (!fas->smoothd) {
1185     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr);
1186   }
1187   if (!fas->smoothu) {
1188     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu);CHKERRQ(ierr);
1189   }
1190   *smooth = fas->smoothu;
1191   PetscFunctionReturn(0);
1192 }
1193 
1194 #undef __FUNCT__
1195 #define __FUNCT__ "SNESFASGetCoarseSolve"
1196 /*@
1197   SNESFASGetCoarseSolve - Gets the coarsest solver.
1198 
1199   Input Parameters:
1200 . snes - the multigrid context
1201 
1202   Output Parameters:
1203 . coarse - the coarse-level solver
1204 
1205   Level: advanced
1206 
1207 .keywords: FAS, MG, get, multigrid, solver, coarse
1208 .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1209 @*/
1210 PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
1211 {
1212   SNES_FAS       *fas;
1213   PetscErrorCode ierr;
1214   SNES           levelsnes;
1215 
1216   PetscFunctionBegin;
1217   ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr);
1218   fas  = (SNES_FAS*)levelsnes->data;
1219   /* if the user chooses to differentiate smoothers, create them both at this point */
1220   if (!fas->smoothd) {
1221     ierr = SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd);CHKERRQ(ierr);
1222   }
1223   *coarse = fas->smoothd;
1224   PetscFunctionReturn(0);
1225 }
1226 
1227 #undef __FUNCT__
1228 #define __FUNCT__ "SNESFASFullSetDownSweep"
1229 /*@
1230    SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS
1231 
1232    Logically Collective on SNES
1233 
1234    Input Parameters:
1235 +  snes - the multigrid context
1236 -  swp - whether to downsweep or not
1237 
1238    Options Database Key:
1239 .  -snes_fas_full_downsweep - Sets number of pre-smoothing steps
1240 
1241    Level: advanced
1242 
1243 .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
1244 
1245 .seealso: SNESFASSetNumberSmoothUp()
1246 @*/
1247 PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp)
1248 {
1249   SNES_FAS       *fas =  (SNES_FAS*)snes->data;
1250   PetscErrorCode ierr = 0;
1251 
1252   PetscFunctionBegin;
1253   fas->full_downsweep = swp;
1254   if (fas->next) {
1255     ierr = SNESFASFullSetDownSweep(fas->next,swp);CHKERRQ(ierr);
1256   }
1257   PetscFunctionReturn(0);
1258 }
1259