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