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