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