xref: /petsc/src/snes/impls/fas/fasfunc.c (revision dbbe0bcd3f3a8fbab5a45420dc06f8387e5764c6)
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 @*/
18ab8d36c9SPeter Brune PetscErrorCode  SNESFASSetType(SNES snes,SNESFASType fastype)
19ab8d36c9SPeter Brune {
20f833ba53SLisandro Dalcin   SNES_FAS       *fas;
2122d28d08SBarry Smith 
22ab8d36c9SPeter Brune   PetscFunctionBegin;
23f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
24ab8d36c9SPeter Brune   PetscValidLogicalCollectiveEnum(snes,fastype,2);
25f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
26ab8d36c9SPeter Brune   fas->fastype = fastype;
271baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetType(fas->next, fastype));
28ab8d36c9SPeter Brune   PetscFunctionReturn(0);
29ab8d36c9SPeter Brune }
30ab8d36c9SPeter Brune 
31ab8d36c9SPeter Brune /*@
32ab8d36c9SPeter Brune SNESFASGetType - Sets the update and correction type used for FAS.
33ab8d36c9SPeter Brune 
34ab8d36c9SPeter Brune Logically Collective
35ab8d36c9SPeter Brune 
36ab8d36c9SPeter Brune Input Parameters:
37ab8d36c9SPeter Brune . snes - FAS context
38ab8d36c9SPeter Brune 
39ab8d36c9SPeter Brune Output Parameters:
40ab8d36c9SPeter Brune . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE
41ab8d36c9SPeter Brune 
42583a1250SSatish Balay Level: intermediate
43583a1250SSatish Balay 
44db781477SPatrick Sanan .seealso: `PCMGSetType()`
45ab8d36c9SPeter Brune @*/
46ab8d36c9SPeter Brune PetscErrorCode  SNESFASGetType(SNES snes,SNESFASType *fastype)
47ab8d36c9SPeter Brune {
48f833ba53SLisandro Dalcin   SNES_FAS *fas;
49ab8d36c9SPeter Brune 
50ab8d36c9SPeter Brune   PetscFunctionBegin;
51f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
52ab8d36c9SPeter Brune   PetscValidPointer(fastype, 2);
53f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
54ab8d36c9SPeter Brune   *fastype = fas->fastype;
55ab8d36c9SPeter Brune   PetscFunctionReturn(0);
56ab8d36c9SPeter Brune }
57ab8d36c9SPeter Brune 
58ab8d36c9SPeter Brune /*@C
59ab8d36c9SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
60ab8d36c9SPeter Brune    Must be called before any other FAS routine.
61ab8d36c9SPeter Brune 
62ab8d36c9SPeter Brune    Input Parameters:
63ab8d36c9SPeter Brune +  snes   - the snes context
64ab8d36c9SPeter Brune .  levels - the number of levels
65ab8d36c9SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
662bf68e3eSBarry Smith             problems on smaller sets of processors.
67ab8d36c9SPeter Brune 
68ab8d36c9SPeter Brune    Level: intermediate
69ab8d36c9SPeter Brune 
70ab8d36c9SPeter Brune    Notes:
71ab8d36c9SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
72ab8d36c9SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
73ab8d36c9SPeter Brune 
74db781477SPatrick Sanan .seealso: `SNESFASGetLevels()`
75ab8d36c9SPeter Brune @*/
7622d28d08SBarry Smith PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms)
7722d28d08SBarry Smith {
78ab8d36c9SPeter Brune   PetscInt       i;
79ab8d36c9SPeter Brune   const char     *optionsprefix;
80ab8d36c9SPeter Brune   char           tprefix[128];
81f833ba53SLisandro Dalcin   SNES_FAS       *fas;
82ab8d36c9SPeter Brune   SNES           prevsnes;
83ab8d36c9SPeter Brune   MPI_Comm       comm;
8422d28d08SBarry Smith 
85ab8d36c9SPeter Brune   PetscFunctionBegin;
86f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
87f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
889566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes,&comm));
89ab8d36c9SPeter Brune   if (levels == fas->levels) {
9022d28d08SBarry Smith     if (!comms) PetscFunctionReturn(0);
91ab8d36c9SPeter Brune   }
92ab8d36c9SPeter Brune   /* user has changed the number of levels; reset */
93*dbbe0bcdSBarry Smith   PetscUseTypeMethod(snes,reset);
94ab8d36c9SPeter Brune   /* destroy any coarser levels if necessary */
959566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&fas->next));
960298fd71SBarry Smith   fas->next     = NULL;
970298fd71SBarry Smith   fas->previous = NULL;
98ab8d36c9SPeter Brune   prevsnes      = snes;
99ab8d36c9SPeter Brune   /* setup the finest level */
1009566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
1019566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1));
102ab8d36c9SPeter Brune   for (i = levels - 1; i >= 0; i--) {
103ab8d36c9SPeter Brune     if (comms) comm = comms[i];
104ab8d36c9SPeter Brune     fas->level  = i;
105ab8d36c9SPeter Brune     fas->levels = levels;
106ab8d36c9SPeter Brune     fas->fine   = snes;
1070298fd71SBarry Smith     fas->next   = NULL;
108ab8d36c9SPeter Brune     if (i > 0) {
1099566063dSJacob Faibussowitsch       PetscCall(SNESCreate(comm, &fas->next));
1109566063dSJacob Faibussowitsch       PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
1119566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_cycle_",(int)fas->level));
1129566063dSJacob Faibussowitsch       PetscCall(SNESAppendOptionsPrefix(fas->next,optionsprefix));
1139566063dSJacob Faibussowitsch       PetscCall(SNESAppendOptionsPrefix(fas->next,tprefix));
1149566063dSJacob Faibussowitsch       PetscCall(SNESSetType(fas->next, SNESFAS));
1159566063dSJacob Faibussowitsch       PetscCall(SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs));
1169566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i));
1179566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1));
1181aa26658SKarl Rupp 
119ab8d36c9SPeter Brune       ((SNES_FAS*)fas->next->data)->previous = prevsnes;
1201aa26658SKarl Rupp 
121ab8d36c9SPeter Brune       prevsnes = fas->next;
122ab8d36c9SPeter Brune       fas      = (SNES_FAS*)prevsnes->data;
123ab8d36c9SPeter Brune     }
124ab8d36c9SPeter Brune   }
125ab8d36c9SPeter Brune   PetscFunctionReturn(0);
126ab8d36c9SPeter Brune }
127ab8d36c9SPeter Brune 
128ab8d36c9SPeter Brune /*@
129ab8d36c9SPeter Brune    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
130ab8d36c9SPeter Brune 
131ab8d36c9SPeter Brune    Input Parameter:
132ab8d36c9SPeter Brune .  snes - the nonlinear solver context
133ab8d36c9SPeter Brune 
134ab8d36c9SPeter Brune    Output parameter:
135ab8d36c9SPeter Brune .  levels - the number of levels
136ab8d36c9SPeter Brune 
137ab8d36c9SPeter Brune    Level: advanced
138ab8d36c9SPeter Brune 
139db781477SPatrick Sanan .seealso: `SNESFASSetLevels()`, `PCMGGetLevels()`
140ab8d36c9SPeter Brune @*/
14122d28d08SBarry Smith PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
14222d28d08SBarry Smith {
143f833ba53SLisandro Dalcin   SNES_FAS *fas;
1445fd66863SKarl Rupp 
145ab8d36c9SPeter Brune   PetscFunctionBegin;
146f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
147064a246eSJacob Faibussowitsch   PetscValidIntPointer(levels,2);
148f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
149ab8d36c9SPeter Brune   *levels = fas->levels;
150ab8d36c9SPeter Brune   PetscFunctionReturn(0);
151ab8d36c9SPeter Brune }
152ab8d36c9SPeter Brune 
153ab8d36c9SPeter Brune /*@
154ab8d36c9SPeter Brune    SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
155ab8d36c9SPeter Brune    level of the FAS hierarchy.
156ab8d36c9SPeter Brune 
157ab8d36c9SPeter Brune    Input Parameters:
158ab8d36c9SPeter Brune +  snes    - the multigrid context
159ab8d36c9SPeter Brune    level   - the level to get
160ab8d36c9SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
161ab8d36c9SPeter Brune 
162ab8d36c9SPeter Brune    Level: advanced
163ab8d36c9SPeter Brune 
164db781477SPatrick Sanan .seealso: `SNESFASSetLevels()`, `SNESFASGetLevels()`
165ab8d36c9SPeter Brune @*/
16622d28d08SBarry Smith PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes)
16722d28d08SBarry Smith {
168f833ba53SLisandro Dalcin   SNES_FAS *fas;
169ab8d36c9SPeter Brune   PetscInt i;
170ab8d36c9SPeter Brune 
171ab8d36c9SPeter Brune   PetscFunctionBegin;
172f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
173f833ba53SLisandro Dalcin   PetscValidPointer(lsnes,3);
174f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
17563a3b9bcSJacob 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);
1760b121fc5SBarry Smith   PetscCheck(fas->level ==  fas->levels - 1,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES");
177ab8d36c9SPeter Brune 
178ab8d36c9SPeter Brune   *lsnes = snes;
179ab8d36c9SPeter Brune   for (i = fas->level; i > level; i--) {
180ab8d36c9SPeter Brune     *lsnes = fas->next;
181ab8d36c9SPeter Brune     fas    = (SNES_FAS*)(*lsnes)->data;
182ab8d36c9SPeter Brune   }
18308401ef6SPierre Jolivet   PetscCheck(fas->level == level,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
184ab8d36c9SPeter Brune   PetscFunctionReturn(0);
185ab8d36c9SPeter Brune }
186ab8d36c9SPeter Brune 
187ab8d36c9SPeter Brune /*@
188ab8d36c9SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
189ab8d36c9SPeter Brune    use on all levels.
190ab8d36c9SPeter Brune 
191ab8d36c9SPeter Brune    Logically Collective on SNES
192ab8d36c9SPeter Brune 
193ab8d36c9SPeter Brune    Input Parameters:
194ab8d36c9SPeter Brune +  snes - the multigrid context
195ab8d36c9SPeter Brune -  n    - the number of smoothing steps
196ab8d36c9SPeter Brune 
197ab8d36c9SPeter Brune    Options Database Key:
198ab8d36c9SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
199ab8d36c9SPeter Brune 
200ab8d36c9SPeter Brune    Level: advanced
201ab8d36c9SPeter Brune 
202db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothDown()`
203ab8d36c9SPeter Brune @*/
20422d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
20522d28d08SBarry Smith {
206f833ba53SLisandro Dalcin   SNES_FAS       *fas;
20722d28d08SBarry Smith 
208ab8d36c9SPeter Brune   PetscFunctionBegin;
209f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
210f833ba53SLisandro Dalcin   fas =  (SNES_FAS*)snes->data;
211ab8d36c9SPeter Brune   fas->max_up_it = n;
212656ede7eSPeter Brune   if (!fas->smoothu && fas->level != 0) {
2139566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
214ab8d36c9SPeter Brune   }
2151baa6e33SBarry Smith   if (fas->smoothu) PetscCall(SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs));
2161baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetNumberSmoothUp(fas->next, n));
217ab8d36c9SPeter Brune   PetscFunctionReturn(0);
218ab8d36c9SPeter Brune }
219ab8d36c9SPeter Brune 
220ab8d36c9SPeter Brune /*@
221ab8d36c9SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
222ab8d36c9SPeter Brune    use on all levels.
223ab8d36c9SPeter Brune 
224ab8d36c9SPeter Brune    Logically Collective on SNES
225ab8d36c9SPeter Brune 
226ab8d36c9SPeter Brune    Input Parameters:
227ab8d36c9SPeter Brune +  snes - the multigrid context
228ab8d36c9SPeter Brune -  n    - the number of smoothing steps
229ab8d36c9SPeter Brune 
230ab8d36c9SPeter Brune    Options Database Key:
231ab8d36c9SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
232ab8d36c9SPeter Brune 
233ab8d36c9SPeter Brune    Level: advanced
234ab8d36c9SPeter Brune 
235db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothUp()`
236ab8d36c9SPeter Brune @*/
23722d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
23822d28d08SBarry Smith {
239f833ba53SLisandro Dalcin   SNES_FAS       *fas;
24022d28d08SBarry Smith 
241ab8d36c9SPeter Brune   PetscFunctionBegin;
242f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
243f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
244ab8d36c9SPeter Brune   if (!fas->smoothd) {
2459566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd));
246ab8d36c9SPeter Brune   }
2479566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs));
2481aa26658SKarl Rupp 
249ab8d36c9SPeter Brune   fas->max_down_it = n;
2501baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetNumberSmoothDown(fas->next, n));
251ab8d36c9SPeter Brune   PetscFunctionReturn(0);
252ab8d36c9SPeter Brune }
253ab8d36c9SPeter Brune 
25487f44e3fSPeter Brune /*@
25587f44e3fSPeter Brune    SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep
25687f44e3fSPeter Brune 
25787f44e3fSPeter Brune    Logically Collective on SNES
25887f44e3fSPeter Brune 
25987f44e3fSPeter Brune    Input Parameters:
26087f44e3fSPeter Brune +  snes - the multigrid context
26187f44e3fSPeter Brune -  n    - the number of smoothing steps
26287f44e3fSPeter Brune 
26387f44e3fSPeter Brune    Options Database Key:
26487f44e3fSPeter Brune .  -snes_fas_continuation - sets continuation to true
26587f44e3fSPeter Brune 
26687f44e3fSPeter Brune    Level: advanced
26787f44e3fSPeter Brune 
26895452b02SPatrick Sanan    Notes:
26995452b02SPatrick Sanan     This sets the prefix on the upsweep smoothers to -fas_continuation
27087f44e3fSPeter Brune 
271db781477SPatrick Sanan .seealso: `SNESFAS`
27287f44e3fSPeter Brune @*/
27387f44e3fSPeter Brune PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation)
27487f44e3fSPeter Brune {
27587f44e3fSPeter Brune   const char     *optionsprefix;
27687f44e3fSPeter Brune   char           tprefix[128];
277f833ba53SLisandro Dalcin   SNES_FAS       *fas;
27887f44e3fSPeter Brune 
27987f44e3fSPeter Brune   PetscFunctionBegin;
280f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
281f833ba53SLisandro Dalcin   fas  = (SNES_FAS*)snes->data;
2829566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
28387f44e3fSPeter Brune   if (!fas->smoothu) {
2849566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
28587f44e3fSPeter Brune   }
2869566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(tprefix,"fas_levels_continuation_",sizeof(tprefix)));
2879566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(fas->smoothu, optionsprefix));
2889566063dSJacob Faibussowitsch   PetscCall(SNESAppendOptionsPrefix(fas->smoothu, tprefix));
2899566063dSJacob Faibussowitsch   PetscCall(SNESSetType(fas->smoothu,SNESNEWTONLS));
2909566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100));
29187f44e3fSPeter Brune   fas->continuation = continuation;
2921baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetContinuation(fas->next,continuation));
29387f44e3fSPeter Brune   PetscFunctionReturn(0);
29487f44e3fSPeter Brune }
29587f44e3fSPeter Brune 
296ab8d36c9SPeter Brune /*@
297ab8d36c9SPeter Brune    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
298ab8d36c9SPeter Brune    complicated cycling.
299ab8d36c9SPeter Brune 
300ab8d36c9SPeter Brune    Logically Collective on SNES
301ab8d36c9SPeter Brune 
302ab8d36c9SPeter Brune    Input Parameters:
303ab8d36c9SPeter Brune +  snes   - the multigrid context
304ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
305ab8d36c9SPeter Brune 
306ab8d36c9SPeter Brune    Options Database Key:
30767b8a455SSatish Balay .  -snes_fas_cycles <1,2> - 1 for V-cycle, 2 for W-cycle
308ab8d36c9SPeter Brune 
309ab8d36c9SPeter Brune    Level: advanced
310ab8d36c9SPeter Brune 
311db781477SPatrick Sanan .seealso: `SNESFASSetCyclesOnLevel()`
312ab8d36c9SPeter Brune @*/
31322d28d08SBarry Smith PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
31422d28d08SBarry Smith {
315f833ba53SLisandro Dalcin   SNES_FAS       *fas;
316ab8d36c9SPeter Brune   PetscBool      isFine;
31722d28d08SBarry Smith 
318ab8d36c9SPeter Brune   PetscFunctionBegin;
319f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
3209566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
321f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
322ab8d36c9SPeter Brune   fas->n_cycles = cycles;
3231aa26658SKarl Rupp   if (!isFine) {
3249566063dSJacob Faibussowitsch     PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
3251aa26658SKarl Rupp   }
3261baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASSetCycles(fas->next, cycles));
327ab8d36c9SPeter Brune   PetscFunctionReturn(0);
328ab8d36c9SPeter Brune }
329ab8d36c9SPeter Brune 
330c8c899caSPeter Brune /*@
331c8c899caSPeter Brune    SNESFASSetMonitor - Sets the method-specific cycle monitoring
332c8c899caSPeter Brune 
333c8c899caSPeter Brune    Logically Collective on SNES
334c8c899caSPeter Brune 
335c8c899caSPeter Brune    Input Parameters:
336c8c899caSPeter Brune +  snes   - the FAS context
337d142ab34SLawrence Mitchell .  vf     - viewer and format structure (may be NULL if flg is FALSE)
338c8c899caSPeter Brune -  flg    - monitor or not
339c8c899caSPeter Brune 
340c8c899caSPeter Brune    Level: advanced
341c8c899caSPeter Brune 
342db781477SPatrick Sanan .seealso: `SNESFASSetCyclesOnLevel()`
343c8c899caSPeter Brune @*/
344d142ab34SLawrence Mitchell PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg)
34522d28d08SBarry Smith {
346f833ba53SLisandro Dalcin   SNES_FAS       *fas;
347c8c899caSPeter Brune   PetscBool      isFine;
348f833ba53SLisandro Dalcin   PetscInt       i, levels;
349c8c899caSPeter Brune   SNES           levelsnes;
35022d28d08SBarry Smith 
351c8c899caSPeter Brune   PetscFunctionBegin;
352f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
3539566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
354f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
355f833ba53SLisandro Dalcin   levels = fas->levels;
356c8c899caSPeter Brune   if (isFine) {
357c8c899caSPeter Brune     for (i = 0; i < levels; i++) {
3589566063dSJacob Faibussowitsch       PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
359c8c899caSPeter Brune       fas  = (SNES_FAS*)levelsnes->data;
360c8c899caSPeter Brune       if (flg) {
361c8c899caSPeter Brune         /* set the monitors for the upsmoother and downsmoother */
3629566063dSJacob Faibussowitsch         PetscCall(SNESMonitorCancel(levelsnes));
363d142ab34SLawrence Mitchell         /* Only register destroy on finest level */
3649566063dSJacob Faibussowitsch         PetscCall(SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL)));
3651aa26658SKarl Rupp       } else if (i != fas->levels - 1) {
366c8c899caSPeter Brune         /* unset the monitors on the coarse levels */
3679566063dSJacob Faibussowitsch         PetscCall(SNESMonitorCancel(levelsnes));
368c8c899caSPeter Brune       }
369c8c899caSPeter Brune     }
370c8c899caSPeter Brune   }
371c8c899caSPeter Brune   PetscFunctionReturn(0);
372c8c899caSPeter Brune }
373c8c899caSPeter Brune 
3740dd27c6cSPeter Brune /*@
3750dd27c6cSPeter Brune    SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels
3760dd27c6cSPeter Brune 
3770dd27c6cSPeter Brune    Logically Collective on SNES
3780dd27c6cSPeter Brune 
3790dd27c6cSPeter Brune    Input Parameters:
3800dd27c6cSPeter Brune +  snes   - the FAS context
3810dd27c6cSPeter Brune -  flg    - monitor or not
3820dd27c6cSPeter Brune 
3830dd27c6cSPeter Brune    Level: advanced
3840dd27c6cSPeter Brune 
385db781477SPatrick Sanan .seealso: `SNESFASSetMonitor()`
3860dd27c6cSPeter Brune @*/
3870dd27c6cSPeter Brune PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
3880dd27c6cSPeter Brune {
389f833ba53SLisandro Dalcin   SNES_FAS       *fas;
3900dd27c6cSPeter Brune   PetscBool      isFine;
391f833ba53SLisandro Dalcin   PetscInt       i, levels;
3920dd27c6cSPeter Brune   SNES           levelsnes;
3930dd27c6cSPeter Brune   char           eventname[128];
3940dd27c6cSPeter Brune 
3950dd27c6cSPeter Brune   PetscFunctionBegin;
396f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
3979566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
398f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
399f833ba53SLisandro Dalcin   levels = fas->levels;
4000dd27c6cSPeter Brune   if (isFine) {
4010dd27c6cSPeter Brune     for (i = 0; i < levels; i++) {
4029566063dSJacob Faibussowitsch       PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
4030dd27c6cSPeter Brune       fas  = (SNES_FAS*)levelsnes->data;
4040dd27c6cSPeter Brune       if (flg) {
4059566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname,sizeof(eventname),"FASSetup  %d",(int)i));
4069566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup));
4079566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname,sizeof(eventname),"FASSmooth %d",(int)i));
4089566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve));
4099566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname,sizeof(eventname),"FASResid  %d",(int)i));
4109566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual));
4119566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname,sizeof(eventname),"FASInterp %d",(int)i));
4129566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict));
4130dd27c6cSPeter Brune       } else {
4140298fd71SBarry Smith         fas->eventsmoothsetup    = 0;
4150298fd71SBarry Smith         fas->eventsmoothsolve    = 0;
4160298fd71SBarry Smith         fas->eventresidual       = 0;
4170298fd71SBarry Smith         fas->eventinterprestrict = 0;
4180dd27c6cSPeter Brune       }
4190dd27c6cSPeter Brune     }
4200dd27c6cSPeter Brune   }
4210dd27c6cSPeter Brune   PetscFunctionReturn(0);
4220dd27c6cSPeter Brune }
4230dd27c6cSPeter Brune 
424ab8d36c9SPeter Brune /*
425ab8d36c9SPeter Brune Creates the default smoother type.
426ab8d36c9SPeter Brune 
42704d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
428ab8d36c9SPeter Brune 
429ab8d36c9SPeter Brune  */
43022d28d08SBarry Smith PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
43122d28d08SBarry Smith {
432ab8d36c9SPeter Brune   SNES_FAS       *fas;
433ab8d36c9SPeter Brune   const char     *optionsprefix;
434ab8d36c9SPeter Brune   char           tprefix[128];
435ab8d36c9SPeter Brune   SNES           nsmooth;
43622d28d08SBarry Smith 
437ab8d36c9SPeter Brune   PetscFunctionBegin;
438f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
439064a246eSJacob Faibussowitsch   PetscValidPointer(smooth,2);
440ab8d36c9SPeter Brune   fas  = (SNES_FAS*)snes->data;
4419566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
442ab8d36c9SPeter Brune   /* create the default smoother */
4439566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth));
444ab8d36c9SPeter Brune   if (fas->level == 0) {
4459566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(tprefix,"fas_coarse_",sizeof(tprefix)));
4469566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
4479566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
4489566063dSJacob Faibussowitsch     PetscCall(SNESSetType(nsmooth, SNESNEWTONLS));
4499566063dSJacob Faibussowitsch     PetscCall(SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs));
450ab8d36c9SPeter Brune   } else {
4519566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_",(int)fas->level));
4529566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
4539566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
4549566063dSJacob Faibussowitsch     PetscCall(SNESSetType(nsmooth, SNESNRICHARDSON));
4559566063dSJacob Faibussowitsch     PetscCall(SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs));
456ab8d36c9SPeter Brune   }
4579566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1));
4589566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth));
4599566063dSJacob Faibussowitsch   PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth));
4609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level));
461ab8d36c9SPeter Brune   *smooth = nsmooth;
462ab8d36c9SPeter Brune   PetscFunctionReturn(0);
463ab8d36c9SPeter Brune }
464ab8d36c9SPeter Brune 
465ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */
466ab8d36c9SPeter Brune 
467ab8d36c9SPeter Brune /*@
468ab8d36c9SPeter Brune    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
469ab8d36c9SPeter Brune 
470ab8d36c9SPeter Brune    Logically Collective on SNES
471ab8d36c9SPeter Brune 
472ab8d36c9SPeter Brune    Input Parameters:
473ab8d36c9SPeter Brune +  snes   - the multigrid context
474ab8d36c9SPeter Brune .  level  - the level to set the number of cycles on
475ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
476ab8d36c9SPeter Brune 
477ab8d36c9SPeter Brune    Level: advanced
478ab8d36c9SPeter Brune 
479db781477SPatrick Sanan .seealso: `SNESFASSetCycles()`
480ab8d36c9SPeter Brune @*/
48122d28d08SBarry Smith PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
48222d28d08SBarry Smith {
483f833ba53SLisandro Dalcin   SNES_FAS       *fas;
48422d28d08SBarry Smith 
485ab8d36c9SPeter Brune   PetscFunctionBegin;
486f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
487f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
488ab8d36c9SPeter Brune   fas->n_cycles = cycles;
4899566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
490ab8d36c9SPeter Brune   PetscFunctionReturn(0);
491ab8d36c9SPeter Brune }
492ab8d36c9SPeter Brune 
493ab8d36c9SPeter Brune /*@
494ab8d36c9SPeter Brune    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
495ab8d36c9SPeter Brune 
496ab8d36c9SPeter Brune    Logically Collective on SNES
497ab8d36c9SPeter Brune 
498ab8d36c9SPeter Brune    Input Parameters:
499ab8d36c9SPeter Brune .  snes   - the multigrid context
500ab8d36c9SPeter Brune 
501ab8d36c9SPeter Brune    Output Parameters:
502ab8d36c9SPeter Brune .  smooth - the smoother
503ab8d36c9SPeter Brune 
504ab8d36c9SPeter Brune    Level: advanced
505ab8d36c9SPeter Brune 
506db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmootherDown()`
507ab8d36c9SPeter Brune @*/
508ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
509ab8d36c9SPeter Brune {
510ab8d36c9SPeter Brune   SNES_FAS *fas;
5115fd66863SKarl Rupp 
512ab8d36c9SPeter Brune   PetscFunctionBegin;
513f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
514f833ba53SLisandro Dalcin   PetscValidPointer(smooth,2);
515ab8d36c9SPeter Brune   fas     = (SNES_FAS*)snes->data;
516ab8d36c9SPeter Brune   *smooth = fas->smoothd;
517ab8d36c9SPeter Brune   PetscFunctionReturn(0);
518ab8d36c9SPeter Brune }
519ab8d36c9SPeter Brune /*@
520ab8d36c9SPeter Brune    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
521ab8d36c9SPeter Brune 
522ab8d36c9SPeter Brune    Logically Collective on SNES
523ab8d36c9SPeter Brune 
524ab8d36c9SPeter Brune    Input Parameters:
525ab8d36c9SPeter Brune .  snes   - the multigrid context
526ab8d36c9SPeter Brune 
527ab8d36c9SPeter Brune    Output Parameters:
528ab8d36c9SPeter Brune .  smoothu - the smoother
529ab8d36c9SPeter Brune 
530ab8d36c9SPeter Brune    Notes:
531ab8d36c9SPeter Brune    Returns the downsmoother if no up smoother is available.  This enables transparent
532ab8d36c9SPeter Brune    default behavior in the process of the solve.
533ab8d36c9SPeter Brune 
534ab8d36c9SPeter Brune    Level: advanced
535ab8d36c9SPeter Brune 
536db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmoother()`, `SNESFASCycleGetSmootherDown()`
537ab8d36c9SPeter Brune @*/
538ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
539ab8d36c9SPeter Brune {
540ab8d36c9SPeter Brune   SNES_FAS *fas;
5415fd66863SKarl Rupp 
542ab8d36c9SPeter Brune   PetscFunctionBegin;
543f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
544f833ba53SLisandro Dalcin   PetscValidPointer(smoothu,2);
545ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
5461aa26658SKarl Rupp   if (!fas->smoothu) *smoothu = fas->smoothd;
5471aa26658SKarl Rupp   else *smoothu = fas->smoothu;
548ab8d36c9SPeter Brune   PetscFunctionReturn(0);
549ab8d36c9SPeter Brune }
550ab8d36c9SPeter Brune 
551ab8d36c9SPeter Brune /*@
552ab8d36c9SPeter Brune    SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
553ab8d36c9SPeter Brune 
554ab8d36c9SPeter Brune    Logically Collective on SNES
555ab8d36c9SPeter Brune 
556ab8d36c9SPeter Brune    Input Parameters:
557ab8d36c9SPeter Brune .  snes   - the multigrid context
558ab8d36c9SPeter Brune 
559ab8d36c9SPeter Brune    Output Parameters:
560ab8d36c9SPeter Brune .  smoothd - the smoother
561ab8d36c9SPeter Brune 
562ab8d36c9SPeter Brune    Level: advanced
563ab8d36c9SPeter Brune 
564db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
565ab8d36c9SPeter Brune @*/
566ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
567ab8d36c9SPeter Brune {
568ab8d36c9SPeter Brune   SNES_FAS *fas;
5695fd66863SKarl Rupp 
570ab8d36c9SPeter Brune   PetscFunctionBegin;
571f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
572f833ba53SLisandro Dalcin   PetscValidPointer(smoothd,2);
573ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
574ab8d36c9SPeter Brune   *smoothd = fas->smoothd;
575ab8d36c9SPeter Brune   PetscFunctionReturn(0);
576ab8d36c9SPeter Brune }
577ab8d36c9SPeter Brune 
578ab8d36c9SPeter Brune /*@
579ab8d36c9SPeter Brune    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level
580ab8d36c9SPeter Brune 
581ab8d36c9SPeter Brune    Logically Collective on SNES
582ab8d36c9SPeter Brune 
583ab8d36c9SPeter Brune    Input Parameters:
584ab8d36c9SPeter Brune .  snes   - the multigrid context
585ab8d36c9SPeter Brune 
586ab8d36c9SPeter Brune    Output Parameters:
587ab8d36c9SPeter Brune .  correction - the coarse correction on this level
588ab8d36c9SPeter Brune 
589ab8d36c9SPeter Brune    Notes:
5900298fd71SBarry Smith    Returns NULL on the coarsest level.
591ab8d36c9SPeter Brune 
592ab8d36c9SPeter Brune    Level: advanced
593ab8d36c9SPeter Brune 
594db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
595ab8d36c9SPeter Brune @*/
596ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
597ab8d36c9SPeter Brune {
598ab8d36c9SPeter Brune   SNES_FAS *fas;
5995fd66863SKarl Rupp 
600ab8d36c9SPeter Brune   PetscFunctionBegin;
601f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
602f833ba53SLisandro Dalcin   PetscValidPointer(correction,2);
603ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
604ab8d36c9SPeter Brune   *correction = fas->next;
605ab8d36c9SPeter Brune   PetscFunctionReturn(0);
606ab8d36c9SPeter Brune }
607ab8d36c9SPeter Brune 
608ab8d36c9SPeter Brune /*@
609ab8d36c9SPeter Brune    SNESFASCycleGetInterpolation - Gets the interpolation on this level
610ab8d36c9SPeter Brune 
611ab8d36c9SPeter Brune    Logically Collective on SNES
612ab8d36c9SPeter Brune 
613ab8d36c9SPeter Brune    Input Parameters:
614ab8d36c9SPeter Brune .  snes   - the multigrid context
615ab8d36c9SPeter Brune 
616ab8d36c9SPeter Brune    Output Parameters:
617ab8d36c9SPeter Brune .  mat    - the interpolation operator on this level
618ab8d36c9SPeter Brune 
619ab8d36c9SPeter Brune    Level: developer
620ab8d36c9SPeter Brune 
621db781477SPatrick Sanan .seealso: `SNESFASCycleGetSmootherUp()`, `SNESFASCycleGetSmoother()`
622ab8d36c9SPeter Brune @*/
623ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
624ab8d36c9SPeter Brune {
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->interpolate;
632ab8d36c9SPeter Brune   PetscFunctionReturn(0);
633ab8d36c9SPeter Brune }
634ab8d36c9SPeter Brune 
635ab8d36c9SPeter Brune /*@
636ab8d36c9SPeter Brune    SNESFASCycleGetRestriction - Gets the restriction 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: `SNESFASGetRestriction()`, `SNESFASCycleGetInterpolation()`
649ab8d36c9SPeter Brune @*/
650ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
651ab8d36c9SPeter Brune {
652ab8d36c9SPeter Brune   SNES_FAS *fas;
6535fd66863SKarl Rupp 
654ab8d36c9SPeter Brune   PetscFunctionBegin;
655f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
656f833ba53SLisandro Dalcin   PetscValidPointer(mat,2);
657ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
658ab8d36c9SPeter Brune   *mat = fas->restrct;
659ab8d36c9SPeter Brune   PetscFunctionReturn(0);
660ab8d36c9SPeter Brune }
661ab8d36c9SPeter Brune 
662ab8d36c9SPeter Brune /*@
663ab8d36c9SPeter Brune    SNESFASCycleGetInjection - Gets the injection on this level
664ab8d36c9SPeter Brune 
665ab8d36c9SPeter Brune    Logically Collective on SNES
666ab8d36c9SPeter Brune 
667ab8d36c9SPeter Brune    Input Parameters:
668ab8d36c9SPeter Brune .  snes   - the multigrid context
669ab8d36c9SPeter Brune 
670ab8d36c9SPeter Brune    Output Parameters:
671ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
672ab8d36c9SPeter Brune 
673ab8d36c9SPeter Brune    Level: developer
674ab8d36c9SPeter Brune 
675db781477SPatrick Sanan .seealso: `SNESFASGetInjection()`, `SNESFASCycleGetRestriction()`
676ab8d36c9SPeter Brune @*/
677ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
678ab8d36c9SPeter Brune {
679ab8d36c9SPeter Brune   SNES_FAS *fas;
6805fd66863SKarl Rupp 
681ab8d36c9SPeter Brune   PetscFunctionBegin;
682f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
683f833ba53SLisandro Dalcin   PetscValidPointer(mat,2);
684ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
685ab8d36c9SPeter Brune   *mat = fas->inject;
686ab8d36c9SPeter Brune   PetscFunctionReturn(0);
687ab8d36c9SPeter Brune }
688ab8d36c9SPeter Brune 
689ab8d36c9SPeter Brune /*@
690ab8d36c9SPeter Brune    SNESFASCycleGetRScale - Gets the injection on this level
691ab8d36c9SPeter Brune 
692ab8d36c9SPeter Brune    Logically Collective on SNES
693ab8d36c9SPeter Brune 
694ab8d36c9SPeter Brune    Input Parameters:
695ab8d36c9SPeter Brune .  snes   - the multigrid context
696ab8d36c9SPeter Brune 
697ab8d36c9SPeter Brune    Output Parameters:
698ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
699ab8d36c9SPeter Brune 
700ab8d36c9SPeter Brune    Level: developer
701ab8d36c9SPeter Brune 
702db781477SPatrick Sanan .seealso: `SNESFASCycleGetRestriction()`, `SNESFASGetRScale()`
703ab8d36c9SPeter Brune @*/
704ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
705ab8d36c9SPeter Brune {
706ab8d36c9SPeter Brune   SNES_FAS *fas;
7075fd66863SKarl Rupp 
708ab8d36c9SPeter Brune   PetscFunctionBegin;
709f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
710f833ba53SLisandro Dalcin   PetscValidPointer(vec,2);
711ab8d36c9SPeter Brune   fas  = (SNES_FAS*)snes->data;
712ab8d36c9SPeter Brune   *vec = fas->rscale;
713ab8d36c9SPeter Brune   PetscFunctionReturn(0);
714ab8d36c9SPeter Brune }
715ab8d36c9SPeter Brune 
716ab8d36c9SPeter Brune /*@
717ab8d36c9SPeter Brune    SNESFASCycleIsFine - Determines if a given cycle is the fine level.
718ab8d36c9SPeter Brune 
719ab8d36c9SPeter Brune    Logically Collective on SNES
720ab8d36c9SPeter Brune 
721ab8d36c9SPeter Brune    Input Parameters:
722ab8d36c9SPeter Brune .  snes   - the FAS context
723ab8d36c9SPeter Brune 
724ab8d36c9SPeter Brune    Output Parameters:
725ab8d36c9SPeter Brune .  flg - indicates if this is the fine level or not
726ab8d36c9SPeter Brune 
727ab8d36c9SPeter Brune    Level: advanced
728ab8d36c9SPeter Brune 
729db781477SPatrick Sanan .seealso: `SNESFASSetLevels()`
730ab8d36c9SPeter Brune @*/
731ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
732ab8d36c9SPeter Brune {
733ab8d36c9SPeter Brune   SNES_FAS *fas;
7345fd66863SKarl Rupp 
735ab8d36c9SPeter Brune   PetscFunctionBegin;
736f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
737534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
738ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
7391aa26658SKarl Rupp   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
7401aa26658SKarl Rupp   else *flg = PETSC_FALSE;
741ab8d36c9SPeter Brune   PetscFunctionReturn(0);
742ab8d36c9SPeter Brune }
743ab8d36c9SPeter Brune 
744ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */
745ab8d36c9SPeter Brune 
746ab8d36c9SPeter Brune /*@
747ab8d36c9SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
748ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
749ab8d36c9SPeter Brune 
750ab8d36c9SPeter Brune    Input Parameters:
751ab8d36c9SPeter Brune +  snes      - the multigrid context
752ab8d36c9SPeter Brune .  mat       - the interpolation operator
753ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
754ab8d36c9SPeter Brune 
755ab8d36c9SPeter Brune    Level: advanced
756ab8d36c9SPeter Brune 
757ab8d36c9SPeter Brune    Notes:
758ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the restriction
759ab8d36c9SPeter Brune     for the same level.
760ab8d36c9SPeter Brune 
761ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
762ab8d36c9SPeter Brune     out from the matrix size which one it is.
763ab8d36c9SPeter Brune 
764db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`, `SNESFASSetRScale()`
765ab8d36c9SPeter Brune @*/
7660adebc6cSBarry Smith PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
7670adebc6cSBarry Smith {
76822d28d08SBarry Smith   SNES_FAS       *fas;
769ab8d36c9SPeter Brune   SNES           levelsnes;
77022d28d08SBarry Smith 
771ab8d36c9SPeter Brune   PetscFunctionBegin;
772f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
773f833ba53SLisandro Dalcin   if (mat) PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
7749566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
775ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
7769566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
7779566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->interpolate));
778ab8d36c9SPeter Brune   fas->interpolate = mat;
779ab8d36c9SPeter Brune   PetscFunctionReturn(0);
780ab8d36c9SPeter Brune }
781ab8d36c9SPeter Brune 
782ab8d36c9SPeter Brune /*@
783ab8d36c9SPeter Brune    SNESFASGetInterpolation - Gets the matrix used to calculate the
784ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
785ab8d36c9SPeter Brune 
786ab8d36c9SPeter Brune    Input Parameters:
787ab8d36c9SPeter Brune +  snes      - the multigrid context
788ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
789ab8d36c9SPeter Brune 
790ab8d36c9SPeter Brune    Output Parameters:
791ab8d36c9SPeter Brune .  mat       - the interpolation operator
792ab8d36c9SPeter Brune 
793ab8d36c9SPeter Brune    Level: advanced
794ab8d36c9SPeter Brune 
795db781477SPatrick Sanan .seealso: `SNESFASSetInterpolation()`, `SNESFASGetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetRScale()`
796ab8d36c9SPeter Brune @*/
79722d28d08SBarry Smith PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
79822d28d08SBarry Smith {
79922d28d08SBarry Smith   SNES_FAS       *fas;
800ab8d36c9SPeter Brune   SNES           levelsnes;
80122d28d08SBarry Smith 
802ab8d36c9SPeter Brune   PetscFunctionBegin;
803f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
804f833ba53SLisandro Dalcin   PetscValidPointer(mat,3);
8059566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
806ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
807ab8d36c9SPeter Brune   *mat = fas->interpolate;
808ab8d36c9SPeter Brune   PetscFunctionReturn(0);
809ab8d36c9SPeter Brune }
810ab8d36c9SPeter Brune 
811ab8d36c9SPeter Brune /*@
812ab8d36c9SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
813ab8d36c9SPeter Brune    from level l to l-1.
814ab8d36c9SPeter Brune 
815ab8d36c9SPeter Brune    Input Parameters:
816ab8d36c9SPeter Brune +  snes  - the multigrid context
817ab8d36c9SPeter Brune .  mat   - the restriction matrix
818ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
819ab8d36c9SPeter Brune 
820ab8d36c9SPeter Brune    Level: advanced
821ab8d36c9SPeter Brune 
822ab8d36c9SPeter Brune    Notes:
823ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the interpolation
824ab8d36c9SPeter Brune     for the same level.
825ab8d36c9SPeter Brune 
826ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
827ab8d36c9SPeter Brune     out from the matrix size which one it is.
828ab8d36c9SPeter Brune 
829ab8d36c9SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
830ab8d36c9SPeter Brune     is used.
831ab8d36c9SPeter Brune 
832db781477SPatrick Sanan .seealso: `SNESFASSetInterpolation()`, `SNESFASSetInjection()`
833ab8d36c9SPeter Brune @*/
83422d28d08SBarry Smith PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
83522d28d08SBarry Smith {
83622d28d08SBarry Smith   SNES_FAS       *fas;
837ab8d36c9SPeter Brune   SNES           levelsnes;
83822d28d08SBarry Smith 
839ab8d36c9SPeter Brune   PetscFunctionBegin;
840f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
841f833ba53SLisandro Dalcin   if (mat) PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
8429566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
843ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
8449566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
8459566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->restrct));
846ab8d36c9SPeter Brune   fas->restrct = mat;
847ab8d36c9SPeter Brune   PetscFunctionReturn(0);
848ab8d36c9SPeter Brune }
849ab8d36c9SPeter Brune 
850ab8d36c9SPeter Brune /*@
851ab8d36c9SPeter Brune    SNESFASGetRestriction - Gets the matrix used to calculate the
852ab8d36c9SPeter Brune    restriction from l to the l-1th level
853ab8d36c9SPeter Brune 
854ab8d36c9SPeter Brune    Input Parameters:
855ab8d36c9SPeter Brune +  snes      - the multigrid context
856ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
857ab8d36c9SPeter Brune 
858ab8d36c9SPeter Brune    Output Parameters:
859ab8d36c9SPeter Brune .  mat       - the interpolation operator
860ab8d36c9SPeter Brune 
861ab8d36c9SPeter Brune    Level: advanced
862ab8d36c9SPeter Brune 
863db781477SPatrick Sanan .seealso: `SNESFASSetRestriction()`, `SNESFASGetInjection()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()`
864ab8d36c9SPeter Brune @*/
86522d28d08SBarry Smith PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
86622d28d08SBarry Smith {
86722d28d08SBarry Smith   SNES_FAS       *fas;
868ab8d36c9SPeter Brune   SNES           levelsnes;
86922d28d08SBarry Smith 
870ab8d36c9SPeter Brune   PetscFunctionBegin;
871f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
872f833ba53SLisandro Dalcin   PetscValidPointer(mat,3);
8739566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
874ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
875ab8d36c9SPeter Brune   *mat = fas->restrct;
876ab8d36c9SPeter Brune   PetscFunctionReturn(0);
877ab8d36c9SPeter Brune }
878ab8d36c9SPeter Brune 
879ab8d36c9SPeter Brune /*@
880ab8d36c9SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
881ab8d36c9SPeter Brune    from level l to l-1.
882ab8d36c9SPeter Brune 
883ab8d36c9SPeter Brune    Input Parameters:
884ab8d36c9SPeter Brune  +  snes  - the multigrid context
885ab8d36c9SPeter Brune .  mat   - the restriction matrix
886ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
887ab8d36c9SPeter Brune 
888ab8d36c9SPeter Brune    Level: advanced
889ab8d36c9SPeter Brune 
890ab8d36c9SPeter Brune    Notes:
891ab8d36c9SPeter Brune          If you do not set this, the restriction and rscale is used to
892ab8d36c9SPeter Brune    project the solution instead.
893ab8d36c9SPeter Brune 
894db781477SPatrick Sanan .seealso: `SNESFASSetInterpolation()`, `SNESFASSetRestriction()`
895ab8d36c9SPeter Brune @*/
89622d28d08SBarry Smith PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
89722d28d08SBarry Smith {
89822d28d08SBarry Smith   SNES_FAS       *fas;
899ab8d36c9SPeter Brune   SNES           levelsnes;
90022d28d08SBarry Smith 
901ab8d36c9SPeter Brune   PetscFunctionBegin;
902f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
903f833ba53SLisandro Dalcin   if (mat) PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
9049566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
905ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
9069566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
9079566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->inject));
9081aa26658SKarl Rupp 
909ab8d36c9SPeter Brune   fas->inject = mat;
910ab8d36c9SPeter Brune   PetscFunctionReturn(0);
911ab8d36c9SPeter Brune }
912ab8d36c9SPeter Brune 
913ab8d36c9SPeter Brune /*@
914ab8d36c9SPeter Brune    SNESFASGetInjection - Gets the matrix used to calculate the
915ab8d36c9SPeter Brune    injection from l-1 to the lth level
916ab8d36c9SPeter Brune 
917ab8d36c9SPeter Brune    Input Parameters:
918ab8d36c9SPeter Brune +  snes      - the multigrid context
919ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
920ab8d36c9SPeter Brune 
921ab8d36c9SPeter Brune    Output Parameters:
922ab8d36c9SPeter Brune .  mat       - the injection operator
923ab8d36c9SPeter Brune 
924ab8d36c9SPeter Brune    Level: advanced
925ab8d36c9SPeter Brune 
926db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASGetRestriction()`, `SNESFASGetInterpolation()`, `SNESFASGetRScale()`
927ab8d36c9SPeter Brune @*/
92822d28d08SBarry Smith PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
92922d28d08SBarry Smith {
93022d28d08SBarry Smith   SNES_FAS       *fas;
931ab8d36c9SPeter Brune   SNES           levelsnes;
93222d28d08SBarry Smith 
933ab8d36c9SPeter Brune   PetscFunctionBegin;
934f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
935f833ba53SLisandro Dalcin   PetscValidPointer(mat,3);
9369566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
937ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
938ab8d36c9SPeter Brune   *mat = fas->inject;
939ab8d36c9SPeter Brune   PetscFunctionReturn(0);
940ab8d36c9SPeter Brune }
941ab8d36c9SPeter Brune 
942ab8d36c9SPeter Brune /*@
943ab8d36c9SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
944ab8d36c9SPeter Brune    operator from level l to l-1.
945ab8d36c9SPeter Brune 
946ab8d36c9SPeter Brune    Input Parameters:
947ab8d36c9SPeter Brune +  snes   - the multigrid context
948ab8d36c9SPeter Brune .  rscale - the restriction scaling
949ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
950ab8d36c9SPeter Brune 
951ab8d36c9SPeter Brune    Level: advanced
952ab8d36c9SPeter Brune 
953ab8d36c9SPeter Brune    Notes:
954ab8d36c9SPeter Brune          This is only used in the case that the injection is not set.
955ab8d36c9SPeter Brune 
956db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
957ab8d36c9SPeter Brune @*/
95822d28d08SBarry Smith PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
95922d28d08SBarry Smith {
960ab8d36c9SPeter Brune   SNES_FAS       *fas;
961ab8d36c9SPeter Brune   SNES           levelsnes;
96222d28d08SBarry Smith 
963ab8d36c9SPeter Brune   PetscFunctionBegin;
964f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
965f833ba53SLisandro Dalcin   if (rscale) PetscValidHeaderSpecific(rscale,VEC_CLASSID,3);
9669566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
967ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
9689566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)rscale));
9699566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fas->rscale));
970ab8d36c9SPeter Brune   fas->rscale = rscale;
971ab8d36c9SPeter Brune   PetscFunctionReturn(0);
972ab8d36c9SPeter Brune }
973ab8d36c9SPeter Brune 
974ab8d36c9SPeter Brune /*@
975ab8d36c9SPeter Brune    SNESFASGetSmoother - Gets the default smoother on a level.
976ab8d36c9SPeter Brune 
977ab8d36c9SPeter Brune    Input Parameters:
978ab8d36c9SPeter Brune +  snes   - the multigrid context
979ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
980ab8d36c9SPeter Brune 
981ab8d36c9SPeter Brune    Output Parameters:
982ab8d36c9SPeter Brune    smooth  - the smoother
983ab8d36c9SPeter Brune 
984ab8d36c9SPeter Brune    Level: advanced
985ab8d36c9SPeter Brune 
986db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
987ab8d36c9SPeter Brune @*/
98822d28d08SBarry Smith PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
98922d28d08SBarry Smith {
990ab8d36c9SPeter Brune   SNES_FAS       *fas;
991ab8d36c9SPeter Brune   SNES           levelsnes;
99222d28d08SBarry Smith 
993ab8d36c9SPeter Brune   PetscFunctionBegin;
994f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
995f833ba53SLisandro Dalcin   PetscValidPointer(smooth,3);
9969566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
997ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
998ab8d36c9SPeter Brune   if (!fas->smoothd) {
9999566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1000ab8d36c9SPeter Brune   }
1001ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1002ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1003ab8d36c9SPeter Brune }
1004ab8d36c9SPeter Brune 
1005ab8d36c9SPeter Brune /*@
1006ab8d36c9SPeter Brune    SNESFASGetSmootherDown - Gets the downsmoother on a level.
1007ab8d36c9SPeter Brune 
1008ab8d36c9SPeter Brune    Input Parameters:
1009ab8d36c9SPeter Brune +  snes   - the multigrid context
1010ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
1011ab8d36c9SPeter Brune 
1012ab8d36c9SPeter Brune    Output Parameters:
1013ab8d36c9SPeter Brune    smooth  - the smoother
1014ab8d36c9SPeter Brune 
1015ab8d36c9SPeter Brune    Level: advanced
1016ab8d36c9SPeter Brune 
1017db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1018ab8d36c9SPeter Brune @*/
101922d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
102022d28d08SBarry Smith {
1021ab8d36c9SPeter Brune   SNES_FAS       *fas;
1022ab8d36c9SPeter Brune   SNES           levelsnes;
102322d28d08SBarry Smith 
1024ab8d36c9SPeter Brune   PetscFunctionBegin;
1025f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
1026f833ba53SLisandro Dalcin   PetscValidPointer(smooth,3);
10279566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1028ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1029ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1030ab8d36c9SPeter Brune   if (!fas->smoothd) {
10319566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1032ab8d36c9SPeter Brune   }
1033ab8d36c9SPeter Brune   if (!fas->smoothu) {
10349566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
1035ab8d36c9SPeter Brune   }
1036ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1037ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1038ab8d36c9SPeter Brune }
1039ab8d36c9SPeter Brune 
1040ab8d36c9SPeter Brune /*@
1041ab8d36c9SPeter Brune    SNESFASGetSmootherUp - Gets the upsmoother on a level.
1042ab8d36c9SPeter Brune 
1043ab8d36c9SPeter Brune    Input Parameters:
1044ab8d36c9SPeter Brune +  snes   - the multigrid context
1045ab8d36c9SPeter Brune -  level  - the level (0 is coarsest)
1046ab8d36c9SPeter Brune 
1047ab8d36c9SPeter Brune    Output Parameters:
1048ab8d36c9SPeter Brune    smooth  - the smoother
1049ab8d36c9SPeter Brune 
1050ab8d36c9SPeter Brune    Level: advanced
1051ab8d36c9SPeter Brune 
1052db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1053ab8d36c9SPeter Brune @*/
105422d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
105522d28d08SBarry Smith {
1056ab8d36c9SPeter Brune   SNES_FAS       *fas;
1057ab8d36c9SPeter Brune   SNES           levelsnes;
105822d28d08SBarry Smith 
1059ab8d36c9SPeter Brune   PetscFunctionBegin;
1060f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
1061f833ba53SLisandro Dalcin   PetscValidPointer(smooth,3);
10629566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1063ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1064ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1065ab8d36c9SPeter Brune   if (!fas->smoothd) {
10669566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1067ab8d36c9SPeter Brune   }
1068ab8d36c9SPeter Brune   if (!fas->smoothu) {
10699566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
1070ab8d36c9SPeter Brune   }
1071ab8d36c9SPeter Brune   *smooth = fas->smoothu;
1072ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1073ab8d36c9SPeter Brune }
1074d6ad1212SPeter Brune 
1075d6ad1212SPeter Brune /*@
1076d6ad1212SPeter Brune   SNESFASGetCoarseSolve - Gets the coarsest solver.
1077d6ad1212SPeter Brune 
1078d6ad1212SPeter Brune   Input Parameters:
1079a3a80b83SMatthew G. Knepley . snes - the multigrid context
1080d6ad1212SPeter Brune 
1081d6ad1212SPeter Brune   Output Parameters:
1082a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver
1083d6ad1212SPeter Brune 
1084d6ad1212SPeter Brune   Level: advanced
1085d6ad1212SPeter Brune 
1086db781477SPatrick Sanan .seealso: `SNESFASSetInjection()`, `SNESFASSetRestriction()`
1087d6ad1212SPeter Brune @*/
1088a3a80b83SMatthew G. Knepley PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
108922d28d08SBarry Smith {
1090d6ad1212SPeter Brune   SNES_FAS       *fas;
1091d6ad1212SPeter Brune   SNES           levelsnes;
109222d28d08SBarry Smith 
1093d6ad1212SPeter Brune   PetscFunctionBegin;
1094f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
1095064a246eSJacob Faibussowitsch   PetscValidPointer(coarse,2);
10969566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, 0, &levelsnes));
1097d6ad1212SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1098d6ad1212SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1099d6ad1212SPeter Brune   if (!fas->smoothd) {
11009566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1101d6ad1212SPeter Brune   }
1102a3a80b83SMatthew G. Knepley   *coarse = fas->smoothd;
1103d6ad1212SPeter Brune   PetscFunctionReturn(0);
1104d6ad1212SPeter Brune }
1105928e959bSPeter Brune 
1106928e959bSPeter Brune /*@
1107928e959bSPeter Brune    SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS
1108928e959bSPeter Brune 
1109928e959bSPeter Brune    Logically Collective on SNES
1110928e959bSPeter Brune 
1111928e959bSPeter Brune    Input Parameters:
1112928e959bSPeter Brune +  snes - the multigrid context
1113928e959bSPeter Brune -  swp - whether to downsweep or not
1114928e959bSPeter Brune 
1115928e959bSPeter Brune    Options Database Key:
1116928e959bSPeter Brune .  -snes_fas_full_downsweep - Sets number of pre-smoothing steps
1117928e959bSPeter Brune 
1118928e959bSPeter Brune    Level: advanced
1119928e959bSPeter Brune 
1120db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothUp()`
1121928e959bSPeter Brune @*/
1122928e959bSPeter Brune PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp)
1123928e959bSPeter Brune {
1124f833ba53SLisandro Dalcin   SNES_FAS       *fas;
1125928e959bSPeter Brune 
1126928e959bSPeter Brune   PetscFunctionBegin;
1127f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
1128f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
1129928e959bSPeter Brune   fas->full_downsweep = swp;
11301baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASFullSetDownSweep(fas->next,swp));
1131928e959bSPeter Brune   PetscFunctionReturn(0);
1132928e959bSPeter Brune }
11336dfbafefSToby Isaac 
11346dfbafefSToby Isaac /*@
11356dfbafefSToby Isaac    SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles
11366dfbafefSToby Isaac 
11376dfbafefSToby Isaac    Logically Collective on SNES
11386dfbafefSToby Isaac 
11396dfbafefSToby Isaac    Input Parameters:
11406dfbafefSToby Isaac +  snes - the multigrid context
11416dfbafefSToby Isaac -  total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
11426dfbafefSToby Isaac 
11436dfbafefSToby Isaac    Options Database Key:
114467b8a455SSatish Balay .  -snes_fas_full_total - Use total restriction and interpolation on the initial down and up sweeps for the full FAS cycle
11456dfbafefSToby Isaac 
11466dfbafefSToby Isaac    Level: advanced
11476dfbafefSToby Isaac 
11486dfbafefSToby Isaac    Note: this option is only significant if the interpolation of a coarse correction (MatInterpolate()) is significantly different from total solution interpolation (DMInterpolateSolution()).
11496dfbafefSToby Isaac 
1150db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`
11516dfbafefSToby Isaac @*/
11526dfbafefSToby Isaac PetscErrorCode SNESFASFullSetTotal(SNES snes,PetscBool total)
11536dfbafefSToby Isaac {
11546dfbafefSToby Isaac   SNES_FAS       *fas;
11556dfbafefSToby Isaac 
11566dfbafefSToby Isaac   PetscFunctionBegin;
11576dfbafefSToby Isaac   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
11586dfbafefSToby Isaac   fas = (SNES_FAS*)snes->data;
11596dfbafefSToby Isaac   fas->full_total = total;
11601baa6e33SBarry Smith   if (fas->next) PetscCall(SNESFASFullSetTotal(fas->next,total));
11616dfbafefSToby Isaac   PetscFunctionReturn(0);
11626dfbafefSToby Isaac }
11636dfbafefSToby Isaac 
11646dfbafefSToby Isaac /*@
11656dfbafefSToby Isaac    SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles
11666dfbafefSToby Isaac 
11676dfbafefSToby Isaac    Logically Collective on SNES
11686dfbafefSToby Isaac 
11696dfbafefSToby Isaac    Input Parameters:
11706dfbafefSToby Isaac .  snes - the multigrid context
11716dfbafefSToby Isaac 
11726dfbafefSToby Isaac    Output:
11736dfbafefSToby Isaac .  total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
11746dfbafefSToby Isaac 
11756dfbafefSToby Isaac    Level: advanced
11766dfbafefSToby Isaac 
1177db781477SPatrick Sanan .seealso: `SNESFASSetNumberSmoothUp()`, `DMInterpolateSolution()`, `SNESFullSetTotal()`
11786dfbafefSToby Isaac @*/
11796dfbafefSToby Isaac PetscErrorCode SNESFASFullGetTotal(SNES snes,PetscBool *total)
11806dfbafefSToby Isaac {
11816dfbafefSToby Isaac   SNES_FAS       *fas;
11826dfbafefSToby Isaac 
11836dfbafefSToby Isaac   PetscFunctionBegin;
11846dfbafefSToby Isaac   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
11856dfbafefSToby Isaac   fas = (SNES_FAS*)snes->data;
11866dfbafefSToby Isaac   *total = fas->full_total;
11876dfbafefSToby Isaac   PetscFunctionReturn(0);
11886dfbafefSToby Isaac }
1189