xref: /petsc/src/snes/impls/fas/fasfunc.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1a7b5fb5fSBarry Smith #include <../src/snes/impls/fas/fasimpls.h> /*I  "petscsnes.h"  I*/
2ab8d36c9SPeter Brune 
3ab8d36c9SPeter Brune /* -------------- functions called on the fine level -------------- */
4ab8d36c9SPeter Brune 
5ab8d36c9SPeter Brune /*@
6ab8d36c9SPeter Brune     SNESFASSetType - Sets the update and correction type used for FAS.
7ab8d36c9SPeter Brune 
8ab8d36c9SPeter Brune    Logically Collective
9ab8d36c9SPeter Brune 
10ab8d36c9SPeter Brune Input Parameters:
11583a1250SSatish Balay + snes  - FAS context
1234d65b3cSPeter Brune - fastype  - SNES_FAS_ADDITIVE, SNES_FAS_MULTIPLICATIVE, SNES_FAS_FULL, or SNES_FAS_KASKADE
13583a1250SSatish Balay 
14583a1250SSatish Balay Level: intermediate
15ab8d36c9SPeter Brune 
16db781477SPatrick Sanan .seealso: `PCMGSetType()`
17ab8d36c9SPeter Brune @*/
189371c9d4SSatish Balay PetscErrorCode SNESFASSetType(SNES snes, SNESFASType fastype) {
19f833ba53SLisandro Dalcin   SNES_FAS *fas;
2022d28d08SBarry Smith 
21ab8d36c9SPeter Brune   PetscFunctionBegin;
22f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
23ab8d36c9SPeter Brune   PetscValidLogicalCollectiveEnum(snes, fastype, 2);
24f833ba53SLisandro Dalcin   fas          = (SNES_FAS *)snes->data;
25ab8d36c9SPeter Brune   fas->fastype = fastype;
261baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetType(fas->next, fastype));
27ab8d36c9SPeter Brune   PetscFunctionReturn(0);
28ab8d36c9SPeter Brune }
29ab8d36c9SPeter Brune 
30ab8d36c9SPeter Brune /*@
31ab8d36c9SPeter Brune SNESFASGetType - Sets the update and correction type used for FAS.
32ab8d36c9SPeter Brune 
33ab8d36c9SPeter Brune Logically Collective
34ab8d36c9SPeter Brune 
35ab8d36c9SPeter Brune Input Parameters:
36ab8d36c9SPeter Brune . snes - FAS context
37ab8d36c9SPeter Brune 
38ab8d36c9SPeter Brune Output Parameters:
39ab8d36c9SPeter Brune . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE
40ab8d36c9SPeter Brune 
41583a1250SSatish Balay Level: intermediate
42583a1250SSatish Balay 
43db781477SPatrick Sanan .seealso: `PCMGSetType()`
44ab8d36c9SPeter Brune @*/
459371c9d4SSatish Balay PetscErrorCode SNESFASGetType(SNES snes, SNESFASType *fastype) {
46f833ba53SLisandro Dalcin   SNES_FAS *fas;
47ab8d36c9SPeter Brune 
48ab8d36c9SPeter Brune   PetscFunctionBegin;
49f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
50ab8d36c9SPeter Brune   PetscValidPointer(fastype, 2);
51f833ba53SLisandro Dalcin   fas      = (SNES_FAS *)snes->data;
52ab8d36c9SPeter Brune   *fastype = fas->fastype;
53ab8d36c9SPeter Brune   PetscFunctionReturn(0);
54ab8d36c9SPeter Brune }
55ab8d36c9SPeter Brune 
56ab8d36c9SPeter Brune /*@C
57ab8d36c9SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
58ab8d36c9SPeter Brune    Must be called before any other FAS routine.
59ab8d36c9SPeter Brune 
60ab8d36c9SPeter Brune    Input Parameters:
61ab8d36c9SPeter Brune +  snes   - the snes context
62ab8d36c9SPeter Brune .  levels - the number of levels
63ab8d36c9SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
642bf68e3eSBarry Smith             problems on smaller sets of processors.
65ab8d36c9SPeter Brune 
66ab8d36c9SPeter Brune    Level: intermediate
67ab8d36c9SPeter Brune 
68ab8d36c9SPeter Brune    Notes:
69ab8d36c9SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
70ab8d36c9SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
71ab8d36c9SPeter Brune 
72db781477SPatrick Sanan .seealso: `SNESFASGetLevels()`
73ab8d36c9SPeter Brune @*/
749371c9d4SSatish Balay PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms) {
75ab8d36c9SPeter Brune   PetscInt    i;
76ab8d36c9SPeter Brune   const char *optionsprefix;
77ab8d36c9SPeter Brune   char        tprefix[128];
78f833ba53SLisandro Dalcin   SNES_FAS   *fas;
79ab8d36c9SPeter Brune   SNES        prevsnes;
80ab8d36c9SPeter Brune   MPI_Comm    comm;
8122d28d08SBarry Smith 
82ab8d36c9SPeter Brune   PetscFunctionBegin;
83f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
84f833ba53SLisandro Dalcin   fas = (SNES_FAS *)snes->data;
859566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
86ab8d36c9SPeter Brune   if (levels == fas->levels) {
8722d28d08SBarry Smith     if (!comms) PetscFunctionReturn(0);
88ab8d36c9SPeter Brune   }
89ab8d36c9SPeter Brune   /* user has changed the number of levels; reset */
90dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes, reset);
91ab8d36c9SPeter Brune   /* destroy any coarser levels if necessary */
929566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&fas->next));
930298fd71SBarry Smith   fas->next     = NULL;
940298fd71SBarry Smith   fas->previous = NULL;
95ab8d36c9SPeter Brune   prevsnes      = snes;
96ab8d36c9SPeter Brune   /* setup the finest level */
979566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
989566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataSetInt((PetscObject)snes, PetscMGLevelId, levels - 1));
99ab8d36c9SPeter Brune   for (i = levels - 1; i >= 0; i--) {
100ab8d36c9SPeter Brune     if (comms) comm = comms[i];
101ab8d36c9SPeter Brune     fas->level  = i;
102ab8d36c9SPeter Brune     fas->levels = levels;
103ab8d36c9SPeter Brune     fas->fine   = snes;
1040298fd71SBarry Smith     fas->next   = NULL;
105ab8d36c9SPeter Brune     if (i > 0) {
1069566063dSJacob Faibussowitsch       PetscCall(SNESCreate(comm, &fas->next));
1079566063dSJacob Faibussowitsch       PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
1089566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%d_cycle_", (int)fas->level));
1099566063dSJacob Faibussowitsch       PetscCall(SNESAppendOptionsPrefix(fas->next, optionsprefix));
1109566063dSJacob Faibussowitsch       PetscCall(SNESAppendOptionsPrefix(fas->next, tprefix));
1119566063dSJacob Faibussowitsch       PetscCall(SNESSetType(fas->next, SNESFAS));
1129566063dSJacob Faibussowitsch       PetscCall(SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs));
1139566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i));
1149566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetInt((PetscObject)fas->next, PetscMGLevelId, i - 1));
1151aa26658SKarl Rupp 
116ab8d36c9SPeter Brune       ((SNES_FAS *)fas->next->data)->previous = prevsnes;
1171aa26658SKarl Rupp 
118ab8d36c9SPeter Brune       prevsnes = fas->next;
119ab8d36c9SPeter Brune       fas      = (SNES_FAS *)prevsnes->data;
120ab8d36c9SPeter Brune     }
121ab8d36c9SPeter Brune   }
122ab8d36c9SPeter Brune   PetscFunctionReturn(0);
123ab8d36c9SPeter Brune }
124ab8d36c9SPeter Brune 
125ab8d36c9SPeter Brune /*@
126ab8d36c9SPeter Brune    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
127ab8d36c9SPeter Brune 
128ab8d36c9SPeter Brune    Input Parameter:
129ab8d36c9SPeter Brune .  snes - the nonlinear solver context
130ab8d36c9SPeter Brune 
131ab8d36c9SPeter Brune    Output parameter:
132ab8d36c9SPeter Brune .  levels - the number of levels
133ab8d36c9SPeter Brune 
134ab8d36c9SPeter Brune    Level: advanced
135ab8d36c9SPeter Brune 
136db781477SPatrick Sanan .seealso: `SNESFASSetLevels()`, `PCMGGetLevels()`
137ab8d36c9SPeter Brune @*/
1389371c9d4SSatish Balay PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels) {
139f833ba53SLisandro Dalcin   SNES_FAS *fas;
1405fd66863SKarl Rupp 
141ab8d36c9SPeter Brune   PetscFunctionBegin;
142f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
143064a246eSJacob Faibussowitsch   PetscValidIntPointer(levels, 2);
144f833ba53SLisandro Dalcin   fas     = (SNES_FAS *)snes->data;
145ab8d36c9SPeter Brune   *levels = fas->levels;
146ab8d36c9SPeter Brune   PetscFunctionReturn(0);
147ab8d36c9SPeter Brune }
148ab8d36c9SPeter Brune 
149ab8d36c9SPeter Brune /*@
150ab8d36c9SPeter Brune    SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
151ab8d36c9SPeter Brune    level of the FAS hierarchy.
152ab8d36c9SPeter Brune 
153ab8d36c9SPeter Brune    Input Parameters:
154ab8d36c9SPeter Brune +  snes    - the multigrid context
155ab8d36c9SPeter Brune    level   - the level to get
156ab8d36c9SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
157ab8d36c9SPeter Brune 
158ab8d36c9SPeter Brune    Level: advanced
159ab8d36c9SPeter Brune 
160db781477SPatrick Sanan .seealso: `SNESFASSetLevels()`, `SNESFASGetLevels()`
161ab8d36c9SPeter Brune @*/
1629371c9d4SSatish Balay PetscErrorCode SNESFASGetCycleSNES(SNES snes, PetscInt level, SNES *lsnes) {
163f833ba53SLisandro Dalcin   SNES_FAS *fas;
164ab8d36c9SPeter Brune   PetscInt  i;
165ab8d36c9SPeter Brune 
166ab8d36c9SPeter Brune   PetscFunctionBegin;
167f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
168f833ba53SLisandro Dalcin   PetscValidPointer(lsnes, 3);
169f833ba53SLisandro Dalcin   fas = (SNES_FAS *)snes->data;
17063a3b9bcSJacob Faibussowitsch   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);
1710b121fc5SBarry Smith   PetscCheck(fas->level == fas->levels - 1, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_OUTOFRANGE, "SNESFASGetCycleSNES may only be called on the finest-level SNES");
172ab8d36c9SPeter Brune 
173ab8d36c9SPeter Brune   *lsnes = snes;
174ab8d36c9SPeter Brune   for (i = fas->level; i > level; i--) {
175ab8d36c9SPeter Brune     *lsnes = fas->next;
176ab8d36c9SPeter Brune     fas    = (SNES_FAS *)(*lsnes)->data;
177ab8d36c9SPeter Brune   }
17808401ef6SPierre Jolivet   PetscCheck(fas->level == level, PetscObjectComm((PetscObject)snes), PETSC_ERR_PLIB, "SNESFAS level hierarchy corrupt");
179ab8d36c9SPeter Brune   PetscFunctionReturn(0);
180ab8d36c9SPeter Brune }
181ab8d36c9SPeter Brune 
182ab8d36c9SPeter Brune /*@
183ab8d36c9SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
184ab8d36c9SPeter Brune    use on all levels.
185ab8d36c9SPeter Brune 
186ab8d36c9SPeter Brune    Logically Collective on SNES
187ab8d36c9SPeter Brune 
188ab8d36c9SPeter Brune    Input Parameters:
189ab8d36c9SPeter Brune +  snes - the multigrid context
190ab8d36c9SPeter Brune -  n    - the number of smoothing steps
191ab8d36c9SPeter Brune 
192ab8d36c9SPeter Brune    Options Database Key:
193ab8d36c9SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
194ab8d36c9SPeter Brune 
195ab8d36c9SPeter Brune    Level: advanced
196ab8d36c9SPeter Brune 
197db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothDown()`
198ab8d36c9SPeter Brune @*/
1999371c9d4SSatish Balay PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n) {
200f833ba53SLisandro Dalcin   SNES_FAS *fas;
20122d28d08SBarry Smith 
202ab8d36c9SPeter Brune   PetscFunctionBegin;
203f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
204f833ba53SLisandro Dalcin   fas            = (SNES_FAS *)snes->data;
205ab8d36c9SPeter Brune   fas->max_up_it = n;
206*48a46eb9SPierre Jolivet   if (!fas->smoothu && fas->level != 0) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
2071baa6e33SBarry Smith   if (fas->smoothu) PetscCall(SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs));
2081baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetNumberSmoothUp(fas->next, n));
209ab8d36c9SPeter Brune   PetscFunctionReturn(0);
210ab8d36c9SPeter Brune }
211ab8d36c9SPeter Brune 
212ab8d36c9SPeter Brune /*@
213ab8d36c9SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
214ab8d36c9SPeter Brune    use on all levels.
215ab8d36c9SPeter Brune 
216ab8d36c9SPeter Brune    Logically Collective on SNES
217ab8d36c9SPeter Brune 
218ab8d36c9SPeter Brune    Input Parameters:
219ab8d36c9SPeter Brune +  snes - the multigrid context
220ab8d36c9SPeter Brune -  n    - the number of smoothing steps
221ab8d36c9SPeter Brune 
222ab8d36c9SPeter Brune    Options Database Key:
223ab8d36c9SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
224ab8d36c9SPeter Brune 
225ab8d36c9SPeter Brune    Level: advanced
226ab8d36c9SPeter Brune 
227db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothUp()`
228ab8d36c9SPeter Brune @*/
2299371c9d4SSatish Balay PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n) {
230f833ba53SLisandro Dalcin   SNES_FAS *fas;
23122d28d08SBarry Smith 
232ab8d36c9SPeter Brune   PetscFunctionBegin;
233f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
234f833ba53SLisandro Dalcin   fas = (SNES_FAS *)snes->data;
235*48a46eb9SPierre Jolivet   if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd));
2369566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs));
2371aa26658SKarl Rupp 
238ab8d36c9SPeter Brune   fas->max_down_it = n;
2391baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetNumberSmoothDown(fas->next, n));
240ab8d36c9SPeter Brune   PetscFunctionReturn(0);
241ab8d36c9SPeter Brune }
242ab8d36c9SPeter Brune 
24387f44e3fSPeter Brune /*@
24487f44e3fSPeter Brune    SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep
24587f44e3fSPeter Brune 
24687f44e3fSPeter Brune    Logically Collective on SNES
24787f44e3fSPeter Brune 
24887f44e3fSPeter Brune    Input Parameters:
24987f44e3fSPeter Brune +  snes - the multigrid context
25087f44e3fSPeter Brune -  n    - the number of smoothing steps
25187f44e3fSPeter Brune 
25287f44e3fSPeter Brune    Options Database Key:
25387f44e3fSPeter Brune .  -snes_fas_continuation - sets continuation to true
25487f44e3fSPeter Brune 
25587f44e3fSPeter Brune    Level: advanced
25687f44e3fSPeter Brune 
25795452b02SPatrick Sanan    Notes:
25895452b02SPatrick Sanan     This sets the prefix on the upsweep smoothers to -fas_continuation
25987f44e3fSPeter Brune 
260db781477SPatrick Sanan .seealso: `SNESFAS`
26187f44e3fSPeter Brune @*/
2629371c9d4SSatish Balay PetscErrorCode SNESFASSetContinuation(SNES snes, PetscBool continuation) {
26387f44e3fSPeter Brune   const char *optionsprefix;
26487f44e3fSPeter Brune   char        tprefix[128];
265f833ba53SLisandro Dalcin   SNES_FAS   *fas;
26687f44e3fSPeter Brune 
26787f44e3fSPeter Brune   PetscFunctionBegin;
268f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
269f833ba53SLisandro Dalcin   fas = (SNES_FAS *)snes->data;
2709566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
271*48a46eb9SPierre Jolivet   if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
2729566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(tprefix, "fas_levels_continuation_", sizeof(tprefix)));
2739566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(fas->smoothu, optionsprefix));
2749566063dSJacob Faibussowitsch   PetscCall(SNESAppendOptionsPrefix(fas->smoothu, tprefix));
2759566063dSJacob Faibussowitsch   PetscCall(SNESSetType(fas->smoothu, SNESNEWTONLS));
2769566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(fas->smoothu, fas->fine->abstol, fas->fine->rtol, fas->fine->stol, 50, 100));
27787f44e3fSPeter Brune   fas->continuation = continuation;
2781baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetContinuation(fas->next, continuation));
27987f44e3fSPeter Brune   PetscFunctionReturn(0);
28087f44e3fSPeter Brune }
28187f44e3fSPeter Brune 
282ab8d36c9SPeter Brune /*@
283ab8d36c9SPeter Brune    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
284ab8d36c9SPeter Brune    complicated cycling.
285ab8d36c9SPeter Brune 
286ab8d36c9SPeter Brune    Logically Collective on SNES
287ab8d36c9SPeter Brune 
288ab8d36c9SPeter Brune    Input Parameters:
289ab8d36c9SPeter Brune +  snes   - the multigrid context
290ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
291ab8d36c9SPeter Brune 
292ab8d36c9SPeter Brune    Options Database Key:
29367b8a455SSatish Balay .  -snes_fas_cycles <1,2> - 1 for V-cycle, 2 for W-cycle
294ab8d36c9SPeter Brune 
295ab8d36c9SPeter Brune    Level: advanced
296ab8d36c9SPeter Brune 
297db781477SPatrick Sanan .seealso: `SNESFASSetCyclesOnLevel()`
298ab8d36c9SPeter Brune @*/
2999371c9d4SSatish Balay PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles) {
300f833ba53SLisandro Dalcin   SNES_FAS *fas;
301ab8d36c9SPeter Brune   PetscBool isFine;
30222d28d08SBarry Smith 
303ab8d36c9SPeter Brune   PetscFunctionBegin;
304f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
3059566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
306f833ba53SLisandro Dalcin   fas           = (SNES_FAS *)snes->data;
307ab8d36c9SPeter Brune   fas->n_cycles = cycles;
308*48a46eb9SPierre Jolivet   if (!isFine) PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
3091baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetCycles(fas->next, cycles));
310ab8d36c9SPeter Brune   PetscFunctionReturn(0);
311ab8d36c9SPeter Brune }
312ab8d36c9SPeter Brune 
313c8c899caSPeter Brune /*@
314c8c899caSPeter Brune    SNESFASSetMonitor - Sets the method-specific cycle monitoring
315c8c899caSPeter Brune 
316c8c899caSPeter Brune    Logically Collective on SNES
317c8c899caSPeter Brune 
318c8c899caSPeter Brune    Input Parameters:
319c8c899caSPeter Brune +  snes   - the FAS context
320d142ab34SLawrence Mitchell .  vf     - viewer and format structure (may be NULL if flg is FALSE)
321c8c899caSPeter Brune -  flg    - monitor or not
322c8c899caSPeter Brune 
323c8c899caSPeter Brune    Level: advanced
324c8c899caSPeter Brune 
325db781477SPatrick Sanan .seealso: `SNESFASSetCyclesOnLevel()`
326c8c899caSPeter Brune @*/
3279371c9d4SSatish Balay PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg) {
328f833ba53SLisandro Dalcin   SNES_FAS *fas;
329c8c899caSPeter Brune   PetscBool isFine;
330f833ba53SLisandro Dalcin   PetscInt  i, levels;
331c8c899caSPeter Brune   SNES      levelsnes;
33222d28d08SBarry Smith 
333c8c899caSPeter Brune   PetscFunctionBegin;
334f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
3359566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
336f833ba53SLisandro Dalcin   fas    = (SNES_FAS *)snes->data;
337f833ba53SLisandro Dalcin   levels = fas->levels;
338c8c899caSPeter Brune   if (isFine) {
339c8c899caSPeter Brune     for (i = 0; i < levels; i++) {
3409566063dSJacob Faibussowitsch       PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
341c8c899caSPeter Brune       fas = (SNES_FAS *)levelsnes->data;
342c8c899caSPeter Brune       if (flg) {
343c8c899caSPeter Brune         /* set the monitors for the upsmoother and downsmoother */
3449566063dSJacob Faibussowitsch         PetscCall(SNESMonitorCancel(levelsnes));
345d142ab34SLawrence Mitchell         /* Only register destroy on finest level */
3469566063dSJacob Faibussowitsch         PetscCall(SNESMonitorSet(levelsnes, (PetscErrorCode(*)(SNES, PetscInt, PetscReal, void *))SNESMonitorDefault, vf, (!i ? (PetscErrorCode(*)(void **))PetscViewerAndFormatDestroy : NULL)));
3471aa26658SKarl Rupp       } else if (i != fas->levels - 1) {
348c8c899caSPeter Brune         /* unset the monitors on the coarse levels */
3499566063dSJacob Faibussowitsch         PetscCall(SNESMonitorCancel(levelsnes));
350c8c899caSPeter Brune       }
351c8c899caSPeter Brune     }
352c8c899caSPeter Brune   }
353c8c899caSPeter Brune   PetscFunctionReturn(0);
354c8c899caSPeter Brune }
355c8c899caSPeter Brune 
3560dd27c6cSPeter Brune /*@
3570dd27c6cSPeter Brune    SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels
3580dd27c6cSPeter Brune 
3590dd27c6cSPeter Brune    Logically Collective on SNES
3600dd27c6cSPeter Brune 
3610dd27c6cSPeter Brune    Input Parameters:
3620dd27c6cSPeter Brune +  snes   - the FAS context
3630dd27c6cSPeter Brune -  flg    - monitor or not
3640dd27c6cSPeter Brune 
3650dd27c6cSPeter Brune    Level: advanced
3660dd27c6cSPeter Brune 
367db781477SPatrick Sanan .seealso: `SNESFASSetMonitor()`
3680dd27c6cSPeter Brune @*/
3699371c9d4SSatish Balay PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg) {
370f833ba53SLisandro Dalcin   SNES_FAS *fas;
3710dd27c6cSPeter Brune   PetscBool isFine;
372f833ba53SLisandro Dalcin   PetscInt  i, levels;
3730dd27c6cSPeter Brune   SNES      levelsnes;
3740dd27c6cSPeter Brune   char      eventname[128];
3750dd27c6cSPeter Brune 
3760dd27c6cSPeter Brune   PetscFunctionBegin;
377f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
3789566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
379f833ba53SLisandro Dalcin   fas    = (SNES_FAS *)snes->data;
380f833ba53SLisandro Dalcin   levels = fas->levels;
3810dd27c6cSPeter Brune   if (isFine) {
3820dd27c6cSPeter Brune     for (i = 0; i < levels; i++) {
3839566063dSJacob Faibussowitsch       PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
3840dd27c6cSPeter Brune       fas = (SNES_FAS *)levelsnes->data;
3850dd27c6cSPeter Brune       if (flg) {
3869566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSetup  %d", (int)i));
3879566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsetup));
3889566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASSmooth %d", (int)i));
3899566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventsmoothsolve));
3909566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASResid  %d", (int)i));
3919566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventresidual));
3929566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname, sizeof(eventname), "FASInterp %d", (int)i));
3939566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname, ((PetscObject)snes)->classid, &fas->eventinterprestrict));
3940dd27c6cSPeter Brune       } else {
3950298fd71SBarry Smith         fas->eventsmoothsetup    = 0;
3960298fd71SBarry Smith         fas->eventsmoothsolve    = 0;
3970298fd71SBarry Smith         fas->eventresidual       = 0;
3980298fd71SBarry Smith         fas->eventinterprestrict = 0;
3990dd27c6cSPeter Brune       }
4000dd27c6cSPeter Brune     }
4010dd27c6cSPeter Brune   }
4020dd27c6cSPeter Brune   PetscFunctionReturn(0);
4030dd27c6cSPeter Brune }
4040dd27c6cSPeter Brune 
405ab8d36c9SPeter Brune /*
406ab8d36c9SPeter Brune Creates the default smoother type.
407ab8d36c9SPeter Brune 
40804d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
409ab8d36c9SPeter Brune 
410ab8d36c9SPeter Brune  */
4119371c9d4SSatish Balay PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth) {
412ab8d36c9SPeter Brune   SNES_FAS   *fas;
413ab8d36c9SPeter Brune   const char *optionsprefix;
414ab8d36c9SPeter Brune   char        tprefix[128];
415ab8d36c9SPeter Brune   SNES        nsmooth;
41622d28d08SBarry Smith 
417ab8d36c9SPeter Brune   PetscFunctionBegin;
418f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
419064a246eSJacob Faibussowitsch   PetscValidPointer(smooth, 2);
420ab8d36c9SPeter Brune   fas = (SNES_FAS *)snes->data;
4219566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
422ab8d36c9SPeter Brune   /* create the default smoother */
4239566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth));
424ab8d36c9SPeter Brune   if (fas->level == 0) {
4259566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(tprefix, "fas_coarse_", sizeof(tprefix)));
4269566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
4279566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
4289566063dSJacob Faibussowitsch     PetscCall(SNESSetType(nsmooth, SNESNEWTONLS));
4299566063dSJacob Faibussowitsch     PetscCall(SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs));
430ab8d36c9SPeter Brune   } else {
4319566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(tprefix, sizeof(tprefix), "fas_levels_%d_", (int)fas->level));
4329566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
4339566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
4349566063dSJacob Faibussowitsch     PetscCall(SNESSetType(nsmooth, SNESNRICHARDSON));
4359566063dSJacob Faibussowitsch     PetscCall(SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs));
436ab8d36c9SPeter Brune   }
4379566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1));
4389566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)snes, (PetscObject)nsmooth));
4399566063dSJacob Faibussowitsch   PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth));
4409566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataSetInt((PetscObject)nsmooth, PetscMGLevelId, fas->level));
441ab8d36c9SPeter Brune   *smooth = nsmooth;
442ab8d36c9SPeter Brune   PetscFunctionReturn(0);
443ab8d36c9SPeter Brune }
444ab8d36c9SPeter Brune 
445ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */
446ab8d36c9SPeter Brune 
447ab8d36c9SPeter Brune /*@
448ab8d36c9SPeter Brune    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
449ab8d36c9SPeter Brune 
450ab8d36c9SPeter Brune    Logically Collective on SNES
451ab8d36c9SPeter Brune 
452ab8d36c9SPeter Brune    Input Parameters:
453ab8d36c9SPeter Brune +  snes   - the multigrid context
454ab8d36c9SPeter Brune .  level  - the level to set the number of cycles on
455ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
456ab8d36c9SPeter Brune 
457ab8d36c9SPeter Brune    Level: advanced
458ab8d36c9SPeter Brune 
459db781477SPatrick Sanan .seealso: `SNESFASSetCycles()`
460ab8d36c9SPeter Brune @*/
4619371c9d4SSatish Balay PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles) {
462f833ba53SLisandro Dalcin   SNES_FAS *fas;
46322d28d08SBarry Smith 
464ab8d36c9SPeter Brune   PetscFunctionBegin;
465f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
466f833ba53SLisandro Dalcin   fas           = (SNES_FAS *)snes->data;
467ab8d36c9SPeter Brune   fas->n_cycles = cycles;
4689566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
469ab8d36c9SPeter Brune   PetscFunctionReturn(0);
470ab8d36c9SPeter Brune }
471ab8d36c9SPeter Brune 
472ab8d36c9SPeter Brune /*@
473ab8d36c9SPeter Brune    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
474ab8d36c9SPeter Brune 
475ab8d36c9SPeter Brune    Logically Collective on SNES
476ab8d36c9SPeter Brune 
477ab8d36c9SPeter Brune    Input Parameters:
478ab8d36c9SPeter Brune .  snes   - the multigrid context
479ab8d36c9SPeter Brune 
480ab8d36c9SPeter Brune    Output Parameters:
481ab8d36c9SPeter Brune .  smooth - the smoother
482ab8d36c9SPeter Brune 
483ab8d36c9SPeter Brune    Level: advanced
484ab8d36c9SPeter Brune 
485db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()`
486ab8d36c9SPeter Brune @*/
4879371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth) {
488ab8d36c9SPeter Brune   SNES_FAS *fas;
4895fd66863SKarl Rupp 
490ab8d36c9SPeter Brune   PetscFunctionBegin;
491f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
492f833ba53SLisandro Dalcin   PetscValidPointer(smooth, 2);
493ab8d36c9SPeter Brune   fas     = (SNES_FAS *)snes->data;
494ab8d36c9SPeter Brune   *smooth = fas->smoothd;
495ab8d36c9SPeter Brune   PetscFunctionReturn(0);
496ab8d36c9SPeter Brune }
497ab8d36c9SPeter Brune /*@
498ab8d36c9SPeter Brune    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
499ab8d36c9SPeter Brune 
500ab8d36c9SPeter Brune    Logically Collective on SNES
501ab8d36c9SPeter Brune 
502ab8d36c9SPeter Brune    Input Parameters:
503ab8d36c9SPeter Brune .  snes   - the multigrid context
504ab8d36c9SPeter Brune 
505ab8d36c9SPeter Brune    Output Parameters:
506ab8d36c9SPeter Brune .  smoothu - the smoother
507ab8d36c9SPeter Brune 
508ab8d36c9SPeter Brune    Notes:
509ab8d36c9SPeter Brune    Returns the downsmoother if no up smoother is available.  This enables transparent
510ab8d36c9SPeter Brune    default behavior in the process of the solve.
511ab8d36c9SPeter Brune 
512ab8d36c9SPeter Brune    Level: advanced
513ab8d36c9SPeter Brune 
514db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherDown()`
515ab8d36c9SPeter Brune @*/
5169371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu) {
517ab8d36c9SPeter Brune   SNES_FAS *fas;
5185fd66863SKarl Rupp 
519ab8d36c9SPeter Brune   PetscFunctionBegin;
520f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
521f833ba53SLisandro Dalcin   PetscValidPointer(smoothu, 2);
522ab8d36c9SPeter Brune   fas = (SNES_FAS *)snes->data;
5231aa26658SKarl Rupp   if (!fas->smoothu) *smoothu = fas->smoothd;
5241aa26658SKarl Rupp   else *smoothu = fas->smoothu;
525ab8d36c9SPeter Brune   PetscFunctionReturn(0);
526ab8d36c9SPeter Brune }
527ab8d36c9SPeter Brune 
528ab8d36c9SPeter Brune /*@
529ab8d36c9SPeter Brune    SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
530ab8d36c9SPeter Brune 
531ab8d36c9SPeter Brune    Logically Collective on SNES
532ab8d36c9SPeter Brune 
533ab8d36c9SPeter Brune    Input Parameters:
534ab8d36c9SPeter Brune .  snes   - the multigrid context
535ab8d36c9SPeter Brune 
536ab8d36c9SPeter Brune    Output Parameters:
537ab8d36c9SPeter Brune .  smoothd - the smoother
538ab8d36c9SPeter Brune 
539ab8d36c9SPeter Brune    Level: advanced
540ab8d36c9SPeter Brune 
541db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
542ab8d36c9SPeter Brune @*/
5439371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd) {
544ab8d36c9SPeter Brune   SNES_FAS *fas;
5455fd66863SKarl Rupp 
546ab8d36c9SPeter Brune   PetscFunctionBegin;
547f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
548f833ba53SLisandro Dalcin   PetscValidPointer(smoothd, 2);
549ab8d36c9SPeter Brune   fas      = (SNES_FAS *)snes->data;
550ab8d36c9SPeter Brune   *smoothd = fas->smoothd;
551ab8d36c9SPeter Brune   PetscFunctionReturn(0);
552ab8d36c9SPeter Brune }
553ab8d36c9SPeter Brune 
554ab8d36c9SPeter Brune /*@
555ab8d36c9SPeter Brune    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level
556ab8d36c9SPeter Brune 
557ab8d36c9SPeter Brune    Logically Collective on SNES
558ab8d36c9SPeter Brune 
559ab8d36c9SPeter Brune    Input Parameters:
560ab8d36c9SPeter Brune .  snes   - the multigrid context
561ab8d36c9SPeter Brune 
562ab8d36c9SPeter Brune    Output Parameters:
563ab8d36c9SPeter Brune .  correction - the coarse correction on this level
564ab8d36c9SPeter Brune 
565ab8d36c9SPeter Brune    Notes:
5660298fd71SBarry Smith    Returns NULL on the coarsest level.
567ab8d36c9SPeter Brune 
568ab8d36c9SPeter Brune    Level: advanced
569ab8d36c9SPeter Brune 
570db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
571ab8d36c9SPeter Brune @*/
5729371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction) {
573ab8d36c9SPeter Brune   SNES_FAS *fas;
5745fd66863SKarl Rupp 
575ab8d36c9SPeter Brune   PetscFunctionBegin;
576f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
577f833ba53SLisandro Dalcin   PetscValidPointer(correction, 2);
578ab8d36c9SPeter Brune   fas         = (SNES_FAS *)snes->data;
579ab8d36c9SPeter Brune   *correction = fas->next;
580ab8d36c9SPeter Brune   PetscFunctionReturn(0);
581ab8d36c9SPeter Brune }
582ab8d36c9SPeter Brune 
583ab8d36c9SPeter Brune /*@
584ab8d36c9SPeter Brune    SNESFASCycleGetInterpolation - Gets the interpolation on this level
585ab8d36c9SPeter Brune 
586ab8d36c9SPeter Brune    Logically Collective on SNES
587ab8d36c9SPeter Brune 
588ab8d36c9SPeter Brune    Input Parameters:
589ab8d36c9SPeter Brune .  snes   - the multigrid context
590ab8d36c9SPeter Brune 
591ab8d36c9SPeter Brune    Output Parameters:
592ab8d36c9SPeter Brune .  mat    - the interpolation operator on this level
593ab8d36c9SPeter Brune 
594ab8d36c9SPeter Brune    Level: developer
595ab8d36c9SPeter Brune 
596db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
597ab8d36c9SPeter Brune @*/
5989371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat) {
599ab8d36c9SPeter Brune   SNES_FAS *fas;
6005fd66863SKarl Rupp 
601ab8d36c9SPeter Brune   PetscFunctionBegin;
602f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
603f833ba53SLisandro Dalcin   PetscValidPointer(mat, 2);
604ab8d36c9SPeter Brune   fas  = (SNES_FAS *)snes->data;
605ab8d36c9SPeter Brune   *mat = fas->interpolate;
606ab8d36c9SPeter Brune   PetscFunctionReturn(0);
607ab8d36c9SPeter Brune }
608ab8d36c9SPeter Brune 
609ab8d36c9SPeter Brune /*@
610ab8d36c9SPeter Brune    SNESFASCycleGetRestriction - Gets the restriction on this level
611ab8d36c9SPeter Brune 
612ab8d36c9SPeter Brune    Logically Collective on SNES
613ab8d36c9SPeter Brune 
614ab8d36c9SPeter Brune    Input Parameters:
615ab8d36c9SPeter Brune .  snes   - the multigrid context
616ab8d36c9SPeter Brune 
617ab8d36c9SPeter Brune    Output Parameters:
618ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
619ab8d36c9SPeter Brune 
620ab8d36c9SPeter Brune    Level: developer
621ab8d36c9SPeter Brune 
622db781477SPatrick Sanan .seealso: `SNESFASGetRestriction()`, `SNESFASCycleGetInterpolation()`
623ab8d36c9SPeter Brune @*/
6249371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat) {
625ab8d36c9SPeter Brune   SNES_FAS *fas;
6265fd66863SKarl Rupp 
627ab8d36c9SPeter Brune   PetscFunctionBegin;
628f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
629f833ba53SLisandro Dalcin   PetscValidPointer(mat, 2);
630ab8d36c9SPeter Brune   fas  = (SNES_FAS *)snes->data;
631ab8d36c9SPeter Brune   *mat = fas->restrct;
632ab8d36c9SPeter Brune   PetscFunctionReturn(0);
633ab8d36c9SPeter Brune }
634ab8d36c9SPeter Brune 
635ab8d36c9SPeter Brune /*@
636ab8d36c9SPeter Brune    SNESFASCycleGetInjection - Gets the injection on this level
637ab8d36c9SPeter Brune 
638ab8d36c9SPeter Brune    Logically Collective on SNES
639ab8d36c9SPeter Brune 
640ab8d36c9SPeter Brune    Input Parameters:
641ab8d36c9SPeter Brune .  snes   - the multigrid context
642ab8d36c9SPeter Brune 
643ab8d36c9SPeter Brune    Output Parameters:
644ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
645ab8d36c9SPeter Brune 
646ab8d36c9SPeter Brune    Level: developer
647ab8d36c9SPeter Brune 
648db781477SPatrick Sanan .seealso: `SNESFASGetInjection()`, `SNESFASCycleGetRestriction()`
649ab8d36c9SPeter Brune @*/
6509371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat) {
651ab8d36c9SPeter Brune   SNES_FAS *fas;
6525fd66863SKarl Rupp 
653ab8d36c9SPeter Brune   PetscFunctionBegin;
654f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
655f833ba53SLisandro Dalcin   PetscValidPointer(mat, 2);
656ab8d36c9SPeter Brune   fas  = (SNES_FAS *)snes->data;
657ab8d36c9SPeter Brune   *mat = fas->inject;
658ab8d36c9SPeter Brune   PetscFunctionReturn(0);
659ab8d36c9SPeter Brune }
660ab8d36c9SPeter Brune 
661ab8d36c9SPeter Brune /*@
662ab8d36c9SPeter Brune    SNESFASCycleGetRScale - Gets the injection on this level
663ab8d36c9SPeter Brune 
664ab8d36c9SPeter Brune    Logically Collective on SNES
665ab8d36c9SPeter Brune 
666ab8d36c9SPeter Brune    Input Parameters:
667ab8d36c9SPeter Brune .  snes   - the multigrid context
668ab8d36c9SPeter Brune 
669ab8d36c9SPeter Brune    Output Parameters:
670ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
671ab8d36c9SPeter Brune 
672ab8d36c9SPeter Brune    Level: developer
673ab8d36c9SPeter Brune 
674db781477SPatrick Sanan .seealso: `SNESFASCycleGetRestriction()`, `SNESFASGetRScale()`
675ab8d36c9SPeter Brune @*/
6769371c9d4SSatish Balay PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec) {
677ab8d36c9SPeter Brune   SNES_FAS *fas;
6785fd66863SKarl Rupp 
679ab8d36c9SPeter Brune   PetscFunctionBegin;
680f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
681f833ba53SLisandro Dalcin   PetscValidPointer(vec, 2);
682ab8d36c9SPeter Brune   fas  = (SNES_FAS *)snes->data;
683ab8d36c9SPeter Brune   *vec = fas->rscale;
684ab8d36c9SPeter Brune   PetscFunctionReturn(0);
685ab8d36c9SPeter Brune }
686ab8d36c9SPeter Brune 
687ab8d36c9SPeter Brune /*@
688ab8d36c9SPeter Brune    SNESFASCycleIsFine - Determines if a given cycle is the fine level.
689ab8d36c9SPeter Brune 
690ab8d36c9SPeter Brune    Logically Collective on SNES
691ab8d36c9SPeter Brune 
692ab8d36c9SPeter Brune    Input Parameters:
693ab8d36c9SPeter Brune .  snes   - the FAS context
694ab8d36c9SPeter Brune 
695ab8d36c9SPeter Brune    Output Parameters:
696ab8d36c9SPeter Brune .  flg - indicates if this is the fine level or not
697ab8d36c9SPeter Brune 
698ab8d36c9SPeter Brune    Level: advanced
699ab8d36c9SPeter Brune 
700db781477SPatrick Sanan .seealso: `SNESFASSetLevels()`
701ab8d36c9SPeter Brune @*/
7029371c9d4SSatish Balay PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg) {
703ab8d36c9SPeter Brune   SNES_FAS *fas;
7045fd66863SKarl Rupp 
705ab8d36c9SPeter Brune   PetscFunctionBegin;
706f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
707534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg, 2);
708ab8d36c9SPeter Brune   fas = (SNES_FAS *)snes->data;
7091aa26658SKarl Rupp   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
7101aa26658SKarl Rupp   else *flg = PETSC_FALSE;
711ab8d36c9SPeter Brune   PetscFunctionReturn(0);
712ab8d36c9SPeter Brune }
713ab8d36c9SPeter Brune 
714ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */
715ab8d36c9SPeter Brune 
716ab8d36c9SPeter Brune /*@
717ab8d36c9SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
718ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
719ab8d36c9SPeter Brune 
720ab8d36c9SPeter Brune    Input Parameters:
721ab8d36c9SPeter Brune +  snes      - the multigrid context
722ab8d36c9SPeter Brune .  mat       - the interpolation operator
723ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
724ab8d36c9SPeter Brune 
725ab8d36c9SPeter Brune    Level: advanced
726ab8d36c9SPeter Brune 
727ab8d36c9SPeter Brune    Notes:
728ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the restriction
729ab8d36c9SPeter Brune     for the same level.
730ab8d36c9SPeter Brune 
731ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
732ab8d36c9SPeter Brune     out from the matrix size which one it is.
733ab8d36c9SPeter Brune 
734db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`, `SNESFASSetRScale()`
735ab8d36c9SPeter Brune @*/
7369371c9d4SSatish Balay PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat) {
73722d28d08SBarry Smith   SNES_FAS *fas;
738ab8d36c9SPeter Brune   SNES      levelsnes;
73922d28d08SBarry Smith 
740ab8d36c9SPeter Brune   PetscFunctionBegin;
741f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
742f833ba53SLisandro Dalcin   if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
7439566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
744ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
7459566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
7469566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->interpolate));
747ab8d36c9SPeter Brune   fas->interpolate = mat;
748ab8d36c9SPeter Brune   PetscFunctionReturn(0);
749ab8d36c9SPeter Brune }
750ab8d36c9SPeter Brune 
751ab8d36c9SPeter Brune /*@
752ab8d36c9SPeter Brune    SNESFASGetInterpolation - Gets the matrix used to calculate the
753ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
754ab8d36c9SPeter Brune 
755ab8d36c9SPeter Brune    Input Parameters:
756ab8d36c9SPeter Brune +  snes      - the multigrid context
757ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
758ab8d36c9SPeter Brune 
759ab8d36c9SPeter Brune    Output Parameters:
760ab8d36c9SPeter Brune .  mat       - the interpolation operator
761ab8d36c9SPeter Brune 
762ab8d36c9SPeter Brune    Level: advanced
763ab8d36c9SPeter Brune 
764db781477SPatrick Sanan .seealso: `SNESFASSetInterpolation()`, `SNESFASGetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetRScale()`
765ab8d36c9SPeter Brune @*/
7669371c9d4SSatish Balay PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat) {
76722d28d08SBarry Smith   SNES_FAS *fas;
768ab8d36c9SPeter Brune   SNES      levelsnes;
76922d28d08SBarry Smith 
770ab8d36c9SPeter Brune   PetscFunctionBegin;
771f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
772f833ba53SLisandro Dalcin   PetscValidPointer(mat, 3);
7739566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
774ab8d36c9SPeter Brune   fas  = (SNES_FAS *)levelsnes->data;
775ab8d36c9SPeter Brune   *mat = fas->interpolate;
776ab8d36c9SPeter Brune   PetscFunctionReturn(0);
777ab8d36c9SPeter Brune }
778ab8d36c9SPeter Brune 
779ab8d36c9SPeter Brune /*@
780ab8d36c9SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
781ab8d36c9SPeter Brune    from level l to l-1.
782ab8d36c9SPeter Brune 
783ab8d36c9SPeter Brune    Input Parameters:
784ab8d36c9SPeter Brune +  snes  - the multigrid context
785ab8d36c9SPeter Brune .  mat   - the restriction matrix
786ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
787ab8d36c9SPeter Brune 
788ab8d36c9SPeter Brune    Level: advanced
789ab8d36c9SPeter Brune 
790ab8d36c9SPeter Brune    Notes:
791ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the interpolation
792ab8d36c9SPeter Brune     for the same level.
793ab8d36c9SPeter Brune 
794ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
795ab8d36c9SPeter Brune     out from the matrix size which one it is.
796ab8d36c9SPeter Brune 
797ab8d36c9SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
798ab8d36c9SPeter Brune     is used.
799ab8d36c9SPeter Brune 
800db781477SPatrick Sanan .seealso: `SNESFASSetInterpolation()`, `SNESFASSetInjection()`
801ab8d36c9SPeter Brune @*/
8029371c9d4SSatish Balay PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat) {
80322d28d08SBarry Smith   SNES_FAS *fas;
804ab8d36c9SPeter Brune   SNES      levelsnes;
80522d28d08SBarry Smith 
806ab8d36c9SPeter Brune   PetscFunctionBegin;
807f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
808f833ba53SLisandro Dalcin   if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
8099566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
810ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
8119566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
8129566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->restrct));
813ab8d36c9SPeter Brune   fas->restrct = mat;
814ab8d36c9SPeter Brune   PetscFunctionReturn(0);
815ab8d36c9SPeter Brune }
816ab8d36c9SPeter Brune 
817ab8d36c9SPeter Brune /*@
818ab8d36c9SPeter Brune    SNESFASGetRestriction - Gets the matrix used to calculate the
819ab8d36c9SPeter Brune    restriction from l to the l-1th level
820ab8d36c9SPeter Brune 
821ab8d36c9SPeter Brune    Input Parameters:
822ab8d36c9SPeter Brune +  snes      - the multigrid context
823ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
824ab8d36c9SPeter Brune 
825ab8d36c9SPeter Brune    Output Parameters:
826ab8d36c9SPeter Brune .  mat       - the interpolation operator
827ab8d36c9SPeter Brune 
828ab8d36c9SPeter Brune    Level: advanced
829ab8d36c9SPeter Brune 
830db781477SPatrick Sanan .seealso: `SNESFASSetRestriction()`, `SNESFASGetInjection()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()`
831ab8d36c9SPeter Brune @*/
8329371c9d4SSatish Balay PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat) {
83322d28d08SBarry Smith   SNES_FAS *fas;
834ab8d36c9SPeter Brune   SNES      levelsnes;
83522d28d08SBarry Smith 
836ab8d36c9SPeter Brune   PetscFunctionBegin;
837f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
838f833ba53SLisandro Dalcin   PetscValidPointer(mat, 3);
8399566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
840ab8d36c9SPeter Brune   fas  = (SNES_FAS *)levelsnes->data;
841ab8d36c9SPeter Brune   *mat = fas->restrct;
842ab8d36c9SPeter Brune   PetscFunctionReturn(0);
843ab8d36c9SPeter Brune }
844ab8d36c9SPeter Brune 
845ab8d36c9SPeter Brune /*@
846ab8d36c9SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
847ab8d36c9SPeter Brune    from level l to l-1.
848ab8d36c9SPeter Brune 
849ab8d36c9SPeter Brune    Input Parameters:
850ab8d36c9SPeter Brune  +  snes  - the multigrid context
851ab8d36c9SPeter Brune .  mat   - the restriction matrix
852ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
853ab8d36c9SPeter Brune 
854ab8d36c9SPeter Brune    Level: advanced
855ab8d36c9SPeter Brune 
856ab8d36c9SPeter Brune    Notes:
857ab8d36c9SPeter Brune          If you do not set this, the restriction and rscale is used to
858ab8d36c9SPeter Brune    project the solution instead.
859ab8d36c9SPeter Brune 
860db781477SPatrick Sanan .seealso: `SNESFASSetInterpolation()`, `SNESFASSetRestriction()`
861ab8d36c9SPeter Brune @*/
8629371c9d4SSatish Balay PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat) {
86322d28d08SBarry Smith   SNES_FAS *fas;
864ab8d36c9SPeter Brune   SNES      levelsnes;
86522d28d08SBarry Smith 
866ab8d36c9SPeter Brune   PetscFunctionBegin;
867f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
868f833ba53SLisandro Dalcin   if (mat) PetscValidHeaderSpecific(mat, MAT_CLASSID, 3);
8699566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
870ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
8719566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
8729566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->inject));
8731aa26658SKarl Rupp 
874ab8d36c9SPeter Brune   fas->inject = mat;
875ab8d36c9SPeter Brune   PetscFunctionReturn(0);
876ab8d36c9SPeter Brune }
877ab8d36c9SPeter Brune 
878ab8d36c9SPeter Brune /*@
879ab8d36c9SPeter Brune    SNESFASGetInjection - Gets the matrix used to calculate the
880ab8d36c9SPeter Brune    injection from l-1 to the lth level
881ab8d36c9SPeter Brune 
882ab8d36c9SPeter Brune    Input Parameters:
883ab8d36c9SPeter Brune +  snes      - the multigrid context
884ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
885ab8d36c9SPeter Brune 
886ab8d36c9SPeter Brune    Output Parameters:
887ab8d36c9SPeter Brune .  mat       - the injection operator
888ab8d36c9SPeter Brune 
889ab8d36c9SPeter Brune    Level: advanced
890ab8d36c9SPeter Brune 
891db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()`
892ab8d36c9SPeter Brune @*/
8939371c9d4SSatish Balay PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat) {
89422d28d08SBarry Smith   SNES_FAS *fas;
895ab8d36c9SPeter Brune   SNES      levelsnes;
89622d28d08SBarry Smith 
897ab8d36c9SPeter Brune   PetscFunctionBegin;
898f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
899f833ba53SLisandro Dalcin   PetscValidPointer(mat, 3);
9009566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
901ab8d36c9SPeter Brune   fas  = (SNES_FAS *)levelsnes->data;
902ab8d36c9SPeter Brune   *mat = fas->inject;
903ab8d36c9SPeter Brune   PetscFunctionReturn(0);
904ab8d36c9SPeter Brune }
905ab8d36c9SPeter Brune 
906ab8d36c9SPeter Brune /*@
907ab8d36c9SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
908ab8d36c9SPeter Brune    operator from level l to l-1.
909ab8d36c9SPeter Brune 
910ab8d36c9SPeter Brune    Input Parameters:
911ab8d36c9SPeter Brune +  snes   - the multigrid context
912ab8d36c9SPeter Brune .  rscale - the restriction scaling
913ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
914ab8d36c9SPeter Brune 
915ab8d36c9SPeter Brune    Level: advanced
916ab8d36c9SPeter Brune 
917ab8d36c9SPeter Brune    Notes:
918ab8d36c9SPeter Brune          This is only used in the case that the injection is not set.
919ab8d36c9SPeter Brune 
920db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
921ab8d36c9SPeter Brune @*/
9229371c9d4SSatish Balay PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale) {
923ab8d36c9SPeter Brune   SNES_FAS *fas;
924ab8d36c9SPeter Brune   SNES      levelsnes;
92522d28d08SBarry Smith 
926ab8d36c9SPeter Brune   PetscFunctionBegin;
927f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
928f833ba53SLisandro Dalcin   if (rscale) PetscValidHeaderSpecific(rscale, VEC_CLASSID, 3);
9299566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
930ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
9319566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)rscale));
9329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fas->rscale));
933ab8d36c9SPeter Brune   fas->rscale = rscale;
934ab8d36c9SPeter Brune   PetscFunctionReturn(0);
935ab8d36c9SPeter Brune }
936ab8d36c9SPeter Brune 
937ab8d36c9SPeter Brune /*@
938ab8d36c9SPeter Brune    SNESFASGetSmoother - Gets the default smoother on a level.
939ab8d36c9SPeter Brune 
940ab8d36c9SPeter Brune    Input Parameters:
941ab8d36c9SPeter Brune +  snes   - the multigrid context
942ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
943ab8d36c9SPeter Brune 
944ab8d36c9SPeter Brune    Output Parameters:
945ab8d36c9SPeter Brune    smooth  - the smoother
946ab8d36c9SPeter Brune 
947ab8d36c9SPeter Brune    Level: advanced
948ab8d36c9SPeter Brune 
949db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
950ab8d36c9SPeter Brune @*/
9519371c9d4SSatish Balay PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth) {
952ab8d36c9SPeter Brune   SNES_FAS *fas;
953ab8d36c9SPeter Brune   SNES      levelsnes;
95422d28d08SBarry Smith 
955ab8d36c9SPeter Brune   PetscFunctionBegin;
956f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
957f833ba53SLisandro Dalcin   PetscValidPointer(smooth, 3);
9589566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
959ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
960*48a46eb9SPierre Jolivet   if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
961ab8d36c9SPeter Brune   *smooth = fas->smoothd;
962ab8d36c9SPeter Brune   PetscFunctionReturn(0);
963ab8d36c9SPeter Brune }
964ab8d36c9SPeter Brune 
965ab8d36c9SPeter Brune /*@
966ab8d36c9SPeter Brune    SNESFASGetSmootherDown - Gets the downsmoother on a level.
967ab8d36c9SPeter Brune 
968ab8d36c9SPeter Brune    Input Parameters:
969ab8d36c9SPeter Brune +  snes   - the multigrid context
970ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
971ab8d36c9SPeter Brune 
972ab8d36c9SPeter Brune    Output Parameters:
973ab8d36c9SPeter Brune    smooth  - the smoother
974ab8d36c9SPeter Brune 
975ab8d36c9SPeter Brune    Level: advanced
976ab8d36c9SPeter Brune 
977db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
978ab8d36c9SPeter Brune @*/
9799371c9d4SSatish Balay PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth) {
980ab8d36c9SPeter Brune   SNES_FAS *fas;
981ab8d36c9SPeter Brune   SNES      levelsnes;
98222d28d08SBarry Smith 
983ab8d36c9SPeter Brune   PetscFunctionBegin;
984f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
985f833ba53SLisandro Dalcin   PetscValidPointer(smooth, 3);
9869566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
987ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
988ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
989*48a46eb9SPierre Jolivet   if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
990*48a46eb9SPierre Jolivet   if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
991ab8d36c9SPeter Brune   *smooth = fas->smoothd;
992ab8d36c9SPeter Brune   PetscFunctionReturn(0);
993ab8d36c9SPeter Brune }
994ab8d36c9SPeter Brune 
995ab8d36c9SPeter Brune /*@
996ab8d36c9SPeter Brune    SNESFASGetSmootherUp - Gets the upsmoother on a level.
997ab8d36c9SPeter Brune 
998ab8d36c9SPeter Brune    Input Parameters:
999ab8d36c9SPeter Brune +  snes   - the multigrid context
1000ab8d36c9SPeter Brune -  level  - the level (0 is coarsest)
1001ab8d36c9SPeter Brune 
1002ab8d36c9SPeter Brune    Output Parameters:
1003ab8d36c9SPeter Brune    smooth  - the smoother
1004ab8d36c9SPeter Brune 
1005ab8d36c9SPeter Brune    Level: advanced
1006ab8d36c9SPeter Brune 
1007db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1008ab8d36c9SPeter Brune @*/
10099371c9d4SSatish Balay PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth) {
1010ab8d36c9SPeter Brune   SNES_FAS *fas;
1011ab8d36c9SPeter Brune   SNES      levelsnes;
101222d28d08SBarry Smith 
1013ab8d36c9SPeter Brune   PetscFunctionBegin;
1014f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
1015f833ba53SLisandro Dalcin   PetscValidPointer(smooth, 3);
10169566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1017ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1018ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1019*48a46eb9SPierre Jolivet   if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1020*48a46eb9SPierre Jolivet   if (!fas->smoothu) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
1021ab8d36c9SPeter Brune   *smooth = fas->smoothu;
1022ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1023ab8d36c9SPeter Brune }
1024d6ad1212SPeter Brune 
1025d6ad1212SPeter Brune /*@
1026d6ad1212SPeter Brune   SNESFASGetCoarseSolve - Gets the coarsest solver.
1027d6ad1212SPeter Brune 
1028d6ad1212SPeter Brune   Input Parameters:
1029a3a80b83SMatthew G. Knepley . snes - the multigrid context
1030d6ad1212SPeter Brune 
1031d6ad1212SPeter Brune   Output Parameters:
1032a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver
1033d6ad1212SPeter Brune 
1034d6ad1212SPeter Brune   Level: advanced
1035d6ad1212SPeter Brune 
1036db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1037d6ad1212SPeter Brune @*/
10389371c9d4SSatish Balay PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse) {
1039d6ad1212SPeter Brune   SNES_FAS *fas;
1040d6ad1212SPeter Brune   SNES      levelsnes;
104122d28d08SBarry Smith 
1042d6ad1212SPeter Brune   PetscFunctionBegin;
1043f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
1044064a246eSJacob Faibussowitsch   PetscValidPointer(coarse, 2);
10459566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, 0, &levelsnes));
1046d6ad1212SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1047d6ad1212SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1048*48a46eb9SPierre Jolivet   if (!fas->smoothd) PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1049a3a80b83SMatthew G. Knepley   *coarse = fas->smoothd;
1050d6ad1212SPeter Brune   PetscFunctionReturn(0);
1051d6ad1212SPeter Brune }
1052928e959bSPeter Brune 
1053928e959bSPeter Brune /*@
1054928e959bSPeter Brune    SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS
1055928e959bSPeter Brune 
1056928e959bSPeter Brune    Logically Collective on SNES
1057928e959bSPeter Brune 
1058928e959bSPeter Brune    Input Parameters:
1059928e959bSPeter Brune +  snes - the multigrid context
1060928e959bSPeter Brune -  swp - whether to downsweep or not
1061928e959bSPeter Brune 
1062928e959bSPeter Brune    Options Database Key:
1063928e959bSPeter Brune .  -snes_fas_full_downsweep - Sets number of pre-smoothing steps
1064928e959bSPeter Brune 
1065928e959bSPeter Brune    Level: advanced
1066928e959bSPeter Brune 
1067db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothUp()`
1068928e959bSPeter Brune @*/
10699371c9d4SSatish Balay PetscErrorCode SNESFASFullSetDownSweep(SNES snes, PetscBool swp) {
1070f833ba53SLisandro Dalcin   SNES_FAS *fas;
1071928e959bSPeter Brune 
1072928e959bSPeter Brune   PetscFunctionBegin;
1073f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
1074f833ba53SLisandro Dalcin   fas                 = (SNES_FAS *)snes->data;
1075928e959bSPeter Brune   fas->full_downsweep = swp;
10761baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASFullSetDownSweep(fas->next, swp));
1077928e959bSPeter Brune   PetscFunctionReturn(0);
1078928e959bSPeter Brune }
10796dfbafefSToby Isaac 
10806dfbafefSToby Isaac /*@
10816dfbafefSToby Isaac    SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles
10826dfbafefSToby Isaac 
10836dfbafefSToby Isaac    Logically Collective on SNES
10846dfbafefSToby Isaac 
10856dfbafefSToby Isaac    Input Parameters:
10866dfbafefSToby Isaac +  snes - the multigrid context
10876dfbafefSToby Isaac -  total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
10886dfbafefSToby Isaac 
10896dfbafefSToby Isaac    Options Database Key:
109067b8a455SSatish Balay .  -snes_fas_full_total - Use total restriction and interpolation on the initial down and up sweeps for the full FAS cycle
10916dfbafefSToby Isaac 
10926dfbafefSToby Isaac    Level: advanced
10936dfbafefSToby Isaac 
10946dfbafefSToby Isaac    Note: this option is only significant if the interpolation of a coarse correction (MatInterpolate()) is significantly different from total solution interpolation (DMInterpolateSolution()).
10956dfbafefSToby Isaac 
1096db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`
10976dfbafefSToby Isaac @*/
10989371c9d4SSatish Balay PetscErrorCode SNESFASFullSetTotal(SNES snes, PetscBool total) {
10996dfbafefSToby Isaac   SNES_FAS *fas;
11006dfbafefSToby Isaac 
11016dfbafefSToby Isaac   PetscFunctionBegin;
11026dfbafefSToby Isaac   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
11036dfbafefSToby Isaac   fas             = (SNES_FAS *)snes->data;
11046dfbafefSToby Isaac   fas->full_total = total;
11051baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASFullSetTotal(fas->next, total));
11066dfbafefSToby Isaac   PetscFunctionReturn(0);
11076dfbafefSToby Isaac }
11086dfbafefSToby Isaac 
11096dfbafefSToby Isaac /*@
11106dfbafefSToby Isaac    SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles
11116dfbafefSToby Isaac 
11126dfbafefSToby Isaac    Logically Collective on SNES
11136dfbafefSToby Isaac 
11146dfbafefSToby Isaac    Input Parameters:
11156dfbafefSToby Isaac .  snes - the multigrid context
11166dfbafefSToby Isaac 
11176dfbafefSToby Isaac    Output:
11186dfbafefSToby Isaac .  total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
11196dfbafefSToby Isaac 
11206dfbafefSToby Isaac    Level: advanced
11216dfbafefSToby Isaac 
1122db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`, `SNESFullSetTotal()`
11236dfbafefSToby Isaac @*/
11249371c9d4SSatish Balay PetscErrorCode SNESFASFullGetTotal(SNES snes, PetscBool *total) {
11256dfbafefSToby Isaac   SNES_FAS *fas;
11266dfbafefSToby Isaac 
11276dfbafefSToby Isaac   PetscFunctionBegin;
11286dfbafefSToby Isaac   PetscValidHeaderSpecificType(snes, SNES_CLASSID, 1, SNESFAS);
11296dfbafefSToby Isaac   fas    = (SNES_FAS *)snes->data;
11306dfbafefSToby Isaac   *total = fas->full_total;
11316dfbafefSToby Isaac   PetscFunctionReturn(0);
11326dfbafefSToby Isaac }
1133