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