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