xref: /petsc/src/snes/impls/fas/fasfunc.c (revision 08401ef684002a709c6d3db98a0c9f54a8bcf1ec)
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 
16ab8d36c9SPeter Brune .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;
2722d28d08SBarry Smith   if (fas->next) {
289566063dSJacob Faibussowitsch     PetscCall(SNESFASSetType(fas->next, fastype));
2922d28d08SBarry Smith   }
30ab8d36c9SPeter Brune   PetscFunctionReturn(0);
31ab8d36c9SPeter Brune }
32ab8d36c9SPeter Brune 
33ab8d36c9SPeter Brune /*@
34ab8d36c9SPeter Brune SNESFASGetType - Sets the update and correction type used for FAS.
35ab8d36c9SPeter Brune 
36ab8d36c9SPeter Brune Logically Collective
37ab8d36c9SPeter Brune 
38ab8d36c9SPeter Brune Input Parameters:
39ab8d36c9SPeter Brune . snes - FAS context
40ab8d36c9SPeter Brune 
41ab8d36c9SPeter Brune Output Parameters:
42ab8d36c9SPeter Brune . fastype - SNES_FAS_ADDITIVE or SNES_FAS_MULTIPLICATIVE
43ab8d36c9SPeter Brune 
44583a1250SSatish Balay Level: intermediate
45583a1250SSatish Balay 
46ab8d36c9SPeter Brune .seealso: PCMGSetType()
47ab8d36c9SPeter Brune @*/
48ab8d36c9SPeter Brune PetscErrorCode  SNESFASGetType(SNES snes,SNESFASType *fastype)
49ab8d36c9SPeter Brune {
50f833ba53SLisandro Dalcin   SNES_FAS *fas;
51ab8d36c9SPeter Brune 
52ab8d36c9SPeter Brune   PetscFunctionBegin;
53f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
54ab8d36c9SPeter Brune   PetscValidPointer(fastype, 2);
55f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
56ab8d36c9SPeter Brune   *fastype = fas->fastype;
57ab8d36c9SPeter Brune   PetscFunctionReturn(0);
58ab8d36c9SPeter Brune }
59ab8d36c9SPeter Brune 
60ab8d36c9SPeter Brune /*@C
61ab8d36c9SPeter Brune    SNESFASSetLevels - Sets the number of levels to use with FAS.
62ab8d36c9SPeter Brune    Must be called before any other FAS routine.
63ab8d36c9SPeter Brune 
64ab8d36c9SPeter Brune    Input Parameters:
65ab8d36c9SPeter Brune +  snes   - the snes context
66ab8d36c9SPeter Brune .  levels - the number of levels
67ab8d36c9SPeter Brune -  comms  - optional communicators for each level; this is to allow solving the coarser
682bf68e3eSBarry Smith             problems on smaller sets of processors.
69ab8d36c9SPeter Brune 
70ab8d36c9SPeter Brune    Level: intermediate
71ab8d36c9SPeter Brune 
72ab8d36c9SPeter Brune    Notes:
73ab8d36c9SPeter Brune      If the number of levels is one then the multigrid uses the -fas_levels prefix
74ab8d36c9SPeter Brune   for setting the level options rather than the -fas_coarse prefix.
75ab8d36c9SPeter Brune 
76ab8d36c9SPeter Brune .seealso: SNESFASGetLevels()
77ab8d36c9SPeter Brune @*/
7822d28d08SBarry Smith PetscErrorCode SNESFASSetLevels(SNES snes, PetscInt levels, MPI_Comm *comms)
7922d28d08SBarry Smith {
80ab8d36c9SPeter Brune   PetscInt       i;
81ab8d36c9SPeter Brune   const char     *optionsprefix;
82ab8d36c9SPeter Brune   char           tprefix[128];
83f833ba53SLisandro Dalcin   SNES_FAS       *fas;
84ab8d36c9SPeter Brune   SNES           prevsnes;
85ab8d36c9SPeter Brune   MPI_Comm       comm;
8622d28d08SBarry Smith 
87ab8d36c9SPeter Brune   PetscFunctionBegin;
88f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
89f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
909566063dSJacob Faibussowitsch   PetscCall(PetscObjectGetComm((PetscObject)snes,&comm));
91ab8d36c9SPeter Brune   if (levels == fas->levels) {
9222d28d08SBarry Smith     if (!comms) PetscFunctionReturn(0);
93ab8d36c9SPeter Brune   }
94ab8d36c9SPeter Brune   /* user has changed the number of levels; reset */
959566063dSJacob Faibussowitsch   PetscCall((*snes->ops->reset)(snes));
96ab8d36c9SPeter Brune   /* destroy any coarser levels if necessary */
979566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&fas->next));
980298fd71SBarry Smith   fas->next     = NULL;
990298fd71SBarry Smith   fas->previous = NULL;
100ab8d36c9SPeter Brune   prevsnes      = snes;
101ab8d36c9SPeter Brune   /* setup the finest level */
1029566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(snes, &optionsprefix));
1039566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataSetInt((PetscObject) snes, PetscMGLevelId, levels-1));
104ab8d36c9SPeter Brune   for (i = levels - 1; i >= 0; i--) {
105ab8d36c9SPeter Brune     if (comms) comm = comms[i];
106ab8d36c9SPeter Brune     fas->level  = i;
107ab8d36c9SPeter Brune     fas->levels = levels;
108ab8d36c9SPeter Brune     fas->fine   = snes;
1090298fd71SBarry Smith     fas->next   = NULL;
110ab8d36c9SPeter Brune     if (i > 0) {
1119566063dSJacob Faibussowitsch       PetscCall(SNESCreate(comm, &fas->next));
1129566063dSJacob Faibussowitsch       PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
1139566063dSJacob Faibussowitsch       PetscCall(PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_cycle_",(int)fas->level));
1149566063dSJacob Faibussowitsch       PetscCall(SNESAppendOptionsPrefix(fas->next,optionsprefix));
1159566063dSJacob Faibussowitsch       PetscCall(SNESAppendOptionsPrefix(fas->next,tprefix));
1169566063dSJacob Faibussowitsch       PetscCall(SNESSetType(fas->next, SNESFAS));
1179566063dSJacob Faibussowitsch       PetscCall(SNESSetTolerances(fas->next, fas->next->abstol, fas->next->rtol, fas->next->stol, fas->n_cycles, fas->next->max_funcs));
1189566063dSJacob Faibussowitsch       PetscCall(PetscObjectIncrementTabLevel((PetscObject)fas->next, (PetscObject)snes, levels - i));
1199566063dSJacob Faibussowitsch       PetscCall(PetscObjectComposedDataSetInt((PetscObject) fas->next, PetscMGLevelId, i-1));
1201aa26658SKarl Rupp 
121ab8d36c9SPeter Brune       ((SNES_FAS*)fas->next->data)->previous = prevsnes;
1221aa26658SKarl Rupp 
123ab8d36c9SPeter Brune       prevsnes = fas->next;
124ab8d36c9SPeter Brune       fas      = (SNES_FAS*)prevsnes->data;
125ab8d36c9SPeter Brune     }
126ab8d36c9SPeter Brune   }
127ab8d36c9SPeter Brune   PetscFunctionReturn(0);
128ab8d36c9SPeter Brune }
129ab8d36c9SPeter Brune 
130ab8d36c9SPeter Brune /*@
131ab8d36c9SPeter Brune    SNESFASGetLevels - Gets the number of levels in a FAS, including fine and coarse grids
132ab8d36c9SPeter Brune 
133ab8d36c9SPeter Brune    Input Parameter:
134ab8d36c9SPeter Brune .  snes - the nonlinear solver context
135ab8d36c9SPeter Brune 
136ab8d36c9SPeter Brune    Output parameter:
137ab8d36c9SPeter Brune .  levels - the number of levels
138ab8d36c9SPeter Brune 
139ab8d36c9SPeter Brune    Level: advanced
140ab8d36c9SPeter Brune 
141ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), PCMGGetLevels()
142ab8d36c9SPeter Brune @*/
14322d28d08SBarry Smith PetscErrorCode SNESFASGetLevels(SNES snes, PetscInt *levels)
14422d28d08SBarry Smith {
145f833ba53SLisandro Dalcin   SNES_FAS *fas;
1465fd66863SKarl Rupp 
147ab8d36c9SPeter Brune   PetscFunctionBegin;
148f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
149064a246eSJacob Faibussowitsch   PetscValidIntPointer(levels,2);
150f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
151ab8d36c9SPeter Brune   *levels = fas->levels;
152ab8d36c9SPeter Brune   PetscFunctionReturn(0);
153ab8d36c9SPeter Brune }
154ab8d36c9SPeter Brune 
155ab8d36c9SPeter Brune /*@
156ab8d36c9SPeter Brune    SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
157ab8d36c9SPeter Brune    level of the FAS hierarchy.
158ab8d36c9SPeter Brune 
159ab8d36c9SPeter Brune    Input Parameters:
160ab8d36c9SPeter Brune +  snes    - the multigrid context
161ab8d36c9SPeter Brune    level   - the level to get
162ab8d36c9SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
163ab8d36c9SPeter Brune 
164ab8d36c9SPeter Brune    Level: advanced
165ab8d36c9SPeter Brune 
166ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
167ab8d36c9SPeter Brune @*/
16822d28d08SBarry Smith PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes)
16922d28d08SBarry Smith {
170f833ba53SLisandro Dalcin   SNES_FAS *fas;
171ab8d36c9SPeter Brune   PetscInt i;
172ab8d36c9SPeter Brune 
173ab8d36c9SPeter Brune   PetscFunctionBegin;
174f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
175f833ba53SLisandro Dalcin   PetscValidPointer(lsnes,3);
176f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
177*08401ef6SPierre Jolivet   PetscCheck(level <= fas->levels-1,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Requested level %D from SNESFAS containing %D levels",level,fas->levels);
1782c71b3e2SJacob Faibussowitsch   PetscCheckFalse(fas->level !=  fas->levels - 1,PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"SNESFASGetCycleSNES may only be called on the finest-level SNES.",level,fas->level);
179ab8d36c9SPeter Brune 
180ab8d36c9SPeter Brune   *lsnes = snes;
181ab8d36c9SPeter Brune   for (i = fas->level; i > level; i--) {
182ab8d36c9SPeter Brune     *lsnes = fas->next;
183ab8d36c9SPeter Brune     fas    = (SNES_FAS*)(*lsnes)->data;
184ab8d36c9SPeter Brune   }
185*08401ef6SPierre Jolivet   PetscCheck(fas->level == level,PetscObjectComm((PetscObject)snes),PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
186ab8d36c9SPeter Brune   PetscFunctionReturn(0);
187ab8d36c9SPeter Brune }
188ab8d36c9SPeter Brune 
189ab8d36c9SPeter Brune /*@
190ab8d36c9SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
191ab8d36c9SPeter Brune    use on all levels.
192ab8d36c9SPeter Brune 
193ab8d36c9SPeter Brune    Logically Collective on SNES
194ab8d36c9SPeter Brune 
195ab8d36c9SPeter Brune    Input Parameters:
196ab8d36c9SPeter Brune +  snes - the multigrid context
197ab8d36c9SPeter Brune -  n    - the number of smoothing steps
198ab8d36c9SPeter Brune 
199ab8d36c9SPeter Brune    Options Database Key:
200ab8d36c9SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
201ab8d36c9SPeter Brune 
202ab8d36c9SPeter Brune    Level: advanced
203ab8d36c9SPeter Brune 
204ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
205ab8d36c9SPeter Brune @*/
20622d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
20722d28d08SBarry Smith {
208f833ba53SLisandro Dalcin   SNES_FAS       *fas;
20922d28d08SBarry Smith 
210ab8d36c9SPeter Brune   PetscFunctionBegin;
211f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
212f833ba53SLisandro Dalcin   fas =  (SNES_FAS*)snes->data;
213ab8d36c9SPeter Brune   fas->max_up_it = n;
214656ede7eSPeter Brune   if (!fas->smoothu && fas->level != 0) {
2159566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
216ab8d36c9SPeter Brune   }
21722d28d08SBarry Smith   if (fas->smoothu) {
2189566063dSJacob Faibussowitsch     PetscCall(SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs));
21922d28d08SBarry Smith   }
220ab8d36c9SPeter Brune   if (fas->next) {
2219566063dSJacob Faibussowitsch     PetscCall(SNESFASSetNumberSmoothUp(fas->next, n));
222ab8d36c9SPeter Brune   }
223ab8d36c9SPeter Brune   PetscFunctionReturn(0);
224ab8d36c9SPeter Brune }
225ab8d36c9SPeter Brune 
226ab8d36c9SPeter Brune /*@
227ab8d36c9SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
228ab8d36c9SPeter Brune    use on all levels.
229ab8d36c9SPeter Brune 
230ab8d36c9SPeter Brune    Logically Collective on SNES
231ab8d36c9SPeter Brune 
232ab8d36c9SPeter Brune    Input Parameters:
233ab8d36c9SPeter Brune +  snes - the multigrid context
234ab8d36c9SPeter Brune -  n    - the number of smoothing steps
235ab8d36c9SPeter Brune 
236ab8d36c9SPeter Brune    Options Database Key:
237ab8d36c9SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
238ab8d36c9SPeter Brune 
239ab8d36c9SPeter Brune    Level: advanced
240ab8d36c9SPeter Brune 
241ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
242ab8d36c9SPeter Brune @*/
24322d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
24422d28d08SBarry Smith {
245f833ba53SLisandro Dalcin   SNES_FAS       *fas;
24622d28d08SBarry Smith 
247ab8d36c9SPeter Brune   PetscFunctionBegin;
248f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
249f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
250ab8d36c9SPeter Brune   if (!fas->smoothd) {
2519566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd));
252ab8d36c9SPeter Brune   }
2539566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs));
2541aa26658SKarl Rupp 
255ab8d36c9SPeter Brune   fas->max_down_it = n;
256ab8d36c9SPeter Brune   if (fas->next) {
2579566063dSJacob Faibussowitsch     PetscCall(SNESFASSetNumberSmoothDown(fas->next, n));
258ab8d36c9SPeter Brune   }
259ab8d36c9SPeter Brune   PetscFunctionReturn(0);
260ab8d36c9SPeter Brune }
261ab8d36c9SPeter Brune 
26287f44e3fSPeter Brune /*@
26387f44e3fSPeter Brune    SNESFASSetContinuation - Sets the FAS cycle to default to exact Newton solves on the upsweep
26487f44e3fSPeter Brune 
26587f44e3fSPeter Brune    Logically Collective on SNES
26687f44e3fSPeter Brune 
26787f44e3fSPeter Brune    Input Parameters:
26887f44e3fSPeter Brune +  snes - the multigrid context
26987f44e3fSPeter Brune -  n    - the number of smoothing steps
27087f44e3fSPeter Brune 
27187f44e3fSPeter Brune    Options Database Key:
27287f44e3fSPeter Brune .  -snes_fas_continuation - sets continuation to true
27387f44e3fSPeter Brune 
27487f44e3fSPeter Brune    Level: advanced
27587f44e3fSPeter Brune 
27695452b02SPatrick Sanan    Notes:
27795452b02SPatrick Sanan     This sets the prefix on the upsweep smoothers to -fas_continuation
27887f44e3fSPeter Brune 
27987f44e3fSPeter Brune .seealso: SNESFAS
28087f44e3fSPeter Brune @*/
28187f44e3fSPeter Brune PetscErrorCode SNESFASSetContinuation(SNES snes,PetscBool continuation)
28287f44e3fSPeter Brune {
28387f44e3fSPeter Brune   const char     *optionsprefix;
28487f44e3fSPeter Brune   char           tprefix[128];
285f833ba53SLisandro Dalcin   SNES_FAS       *fas;
28687f44e3fSPeter Brune 
28787f44e3fSPeter Brune   PetscFunctionBegin;
288f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
289f833ba53SLisandro Dalcin   fas  = (SNES_FAS*)snes->data;
2909566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
29187f44e3fSPeter Brune   if (!fas->smoothu) {
2929566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu));
29387f44e3fSPeter Brune   }
2949566063dSJacob Faibussowitsch   PetscCall(PetscStrncpy(tprefix,"fas_levels_continuation_",sizeof(tprefix)));
2959566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(fas->smoothu, optionsprefix));
2969566063dSJacob Faibussowitsch   PetscCall(SNESAppendOptionsPrefix(fas->smoothu, tprefix));
2979566063dSJacob Faibussowitsch   PetscCall(SNESSetType(fas->smoothu,SNESNEWTONLS));
2989566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(fas->smoothu,fas->fine->abstol,fas->fine->rtol,fas->fine->stol,50,100));
29987f44e3fSPeter Brune   fas->continuation = continuation;
30087f44e3fSPeter Brune   if (fas->next) {
3019566063dSJacob Faibussowitsch     PetscCall(SNESFASSetContinuation(fas->next,continuation));
30287f44e3fSPeter Brune   }
30387f44e3fSPeter Brune   PetscFunctionReturn(0);
30487f44e3fSPeter Brune }
30587f44e3fSPeter Brune 
306ab8d36c9SPeter Brune /*@
307ab8d36c9SPeter Brune    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
308ab8d36c9SPeter Brune    complicated cycling.
309ab8d36c9SPeter Brune 
310ab8d36c9SPeter Brune    Logically Collective on SNES
311ab8d36c9SPeter Brune 
312ab8d36c9SPeter Brune    Input Parameters:
313ab8d36c9SPeter Brune +  snes   - the multigrid context
314ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
315ab8d36c9SPeter Brune 
316ab8d36c9SPeter Brune    Options Database Key:
31767b8a455SSatish Balay .  -snes_fas_cycles <1,2> - 1 for V-cycle, 2 for W-cycle
318ab8d36c9SPeter Brune 
319ab8d36c9SPeter Brune    Level: advanced
320ab8d36c9SPeter Brune 
321ab8d36c9SPeter Brune .seealso: SNESFASSetCyclesOnLevel()
322ab8d36c9SPeter Brune @*/
32322d28d08SBarry Smith PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
32422d28d08SBarry Smith {
325f833ba53SLisandro Dalcin   SNES_FAS       *fas;
326ab8d36c9SPeter Brune   PetscBool      isFine;
32722d28d08SBarry Smith 
328ab8d36c9SPeter Brune   PetscFunctionBegin;
329f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
3309566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
331f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
332ab8d36c9SPeter Brune   fas->n_cycles = cycles;
3331aa26658SKarl Rupp   if (!isFine) {
3349566063dSJacob Faibussowitsch     PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
3351aa26658SKarl Rupp   }
336ab8d36c9SPeter Brune   if (fas->next) {
3379566063dSJacob Faibussowitsch     PetscCall(SNESFASSetCycles(fas->next, cycles));
338ab8d36c9SPeter Brune   }
339ab8d36c9SPeter Brune   PetscFunctionReturn(0);
340ab8d36c9SPeter Brune }
341ab8d36c9SPeter Brune 
342c8c899caSPeter Brune /*@
343c8c899caSPeter Brune    SNESFASSetMonitor - Sets the method-specific cycle monitoring
344c8c899caSPeter Brune 
345c8c899caSPeter Brune    Logically Collective on SNES
346c8c899caSPeter Brune 
347c8c899caSPeter Brune    Input Parameters:
348c8c899caSPeter Brune +  snes   - the FAS context
349d142ab34SLawrence Mitchell .  vf     - viewer and format structure (may be NULL if flg is FALSE)
350c8c899caSPeter Brune -  flg    - monitor or not
351c8c899caSPeter Brune 
352c8c899caSPeter Brune    Level: advanced
353c8c899caSPeter Brune 
354c8c899caSPeter Brune .seealso: SNESFASSetCyclesOnLevel()
355c8c899caSPeter Brune @*/
356d142ab34SLawrence Mitchell PetscErrorCode SNESFASSetMonitor(SNES snes, PetscViewerAndFormat *vf, PetscBool flg)
35722d28d08SBarry Smith {
358f833ba53SLisandro Dalcin   SNES_FAS       *fas;
359c8c899caSPeter Brune   PetscBool      isFine;
360f833ba53SLisandro Dalcin   PetscInt       i, levels;
361c8c899caSPeter Brune   SNES           levelsnes;
36222d28d08SBarry Smith 
363c8c899caSPeter Brune   PetscFunctionBegin;
364f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
3659566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
366f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
367f833ba53SLisandro Dalcin   levels = fas->levels;
368c8c899caSPeter Brune   if (isFine) {
369c8c899caSPeter Brune     for (i = 0; i < levels; i++) {
3709566063dSJacob Faibussowitsch       PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
371c8c899caSPeter Brune       fas  = (SNES_FAS*)levelsnes->data;
372c8c899caSPeter Brune       if (flg) {
373c8c899caSPeter Brune         /* set the monitors for the upsmoother and downsmoother */
3749566063dSJacob Faibussowitsch         PetscCall(SNESMonitorCancel(levelsnes));
375d142ab34SLawrence Mitchell         /* Only register destroy on finest level */
3769566063dSJacob Faibussowitsch         PetscCall(SNESMonitorSet(levelsnes,(PetscErrorCode (*)(SNES,PetscInt,PetscReal,void*))SNESMonitorDefault,vf,(!i ? (PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy : NULL)));
3771aa26658SKarl Rupp       } else if (i != fas->levels - 1) {
378c8c899caSPeter Brune         /* unset the monitors on the coarse levels */
3799566063dSJacob Faibussowitsch         PetscCall(SNESMonitorCancel(levelsnes));
380c8c899caSPeter Brune       }
381c8c899caSPeter Brune     }
382c8c899caSPeter Brune   }
383c8c899caSPeter Brune   PetscFunctionReturn(0);
384c8c899caSPeter Brune }
385c8c899caSPeter Brune 
3860dd27c6cSPeter Brune /*@
3870dd27c6cSPeter Brune    SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels
3880dd27c6cSPeter Brune 
3890dd27c6cSPeter Brune    Logically Collective on SNES
3900dd27c6cSPeter Brune 
3910dd27c6cSPeter Brune    Input Parameters:
3920dd27c6cSPeter Brune +  snes   - the FAS context
3930dd27c6cSPeter Brune -  flg    - monitor or not
3940dd27c6cSPeter Brune 
3950dd27c6cSPeter Brune    Level: advanced
3960dd27c6cSPeter Brune 
3970dd27c6cSPeter Brune .seealso: SNESFASSetMonitor()
3980dd27c6cSPeter Brune @*/
3990dd27c6cSPeter Brune PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
4000dd27c6cSPeter Brune {
401f833ba53SLisandro Dalcin   SNES_FAS       *fas;
4020dd27c6cSPeter Brune   PetscBool      isFine;
403f833ba53SLisandro Dalcin   PetscInt       i, levels;
4040dd27c6cSPeter Brune   SNES           levelsnes;
4050dd27c6cSPeter Brune   char           eventname[128];
4060dd27c6cSPeter Brune 
4070dd27c6cSPeter Brune   PetscFunctionBegin;
408f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
4099566063dSJacob Faibussowitsch   PetscCall(SNESFASCycleIsFine(snes, &isFine));
410f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
411f833ba53SLisandro Dalcin   levels = fas->levels;
4120dd27c6cSPeter Brune   if (isFine) {
4130dd27c6cSPeter Brune     for (i = 0; i < levels; i++) {
4149566063dSJacob Faibussowitsch       PetscCall(SNESFASGetCycleSNES(snes, i, &levelsnes));
4150dd27c6cSPeter Brune       fas  = (SNES_FAS*)levelsnes->data;
4160dd27c6cSPeter Brune       if (flg) {
4179566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname,sizeof(eventname),"FASSetup  %d",(int)i));
4189566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup));
4199566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname,sizeof(eventname),"FASSmooth %d",(int)i));
4209566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve));
4219566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname,sizeof(eventname),"FASResid  %d",(int)i));
4229566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual));
4239566063dSJacob Faibussowitsch         PetscCall(PetscSNPrintf(eventname,sizeof(eventname),"FASInterp %d",(int)i));
4249566063dSJacob Faibussowitsch         PetscCall(PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict));
4250dd27c6cSPeter Brune       } else {
4260298fd71SBarry Smith         fas->eventsmoothsetup    = 0;
4270298fd71SBarry Smith         fas->eventsmoothsolve    = 0;
4280298fd71SBarry Smith         fas->eventresidual       = 0;
4290298fd71SBarry Smith         fas->eventinterprestrict = 0;
4300dd27c6cSPeter Brune       }
4310dd27c6cSPeter Brune     }
4320dd27c6cSPeter Brune   }
4330dd27c6cSPeter Brune   PetscFunctionReturn(0);
4340dd27c6cSPeter Brune }
4350dd27c6cSPeter Brune 
436ab8d36c9SPeter Brune /*
437ab8d36c9SPeter Brune Creates the default smoother type.
438ab8d36c9SPeter Brune 
43904d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
440ab8d36c9SPeter Brune 
441ab8d36c9SPeter Brune  */
44222d28d08SBarry Smith PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
44322d28d08SBarry Smith {
444ab8d36c9SPeter Brune   SNES_FAS       *fas;
445ab8d36c9SPeter Brune   const char     *optionsprefix;
446ab8d36c9SPeter Brune   char           tprefix[128];
447ab8d36c9SPeter Brune   SNES           nsmooth;
44822d28d08SBarry Smith 
449ab8d36c9SPeter Brune   PetscFunctionBegin;
450f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
451064a246eSJacob Faibussowitsch   PetscValidPointer(smooth,2);
452ab8d36c9SPeter Brune   fas  = (SNES_FAS*)snes->data;
4539566063dSJacob Faibussowitsch   PetscCall(SNESGetOptionsPrefix(fas->fine, &optionsprefix));
454ab8d36c9SPeter Brune   /* create the default smoother */
4559566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PetscObjectComm((PetscObject)snes), &nsmooth));
456ab8d36c9SPeter Brune   if (fas->level == 0) {
4579566063dSJacob Faibussowitsch     PetscCall(PetscStrncpy(tprefix,"fas_coarse_",sizeof(tprefix)));
4589566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
4599566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
4609566063dSJacob Faibussowitsch     PetscCall(SNESSetType(nsmooth, SNESNEWTONLS));
4619566063dSJacob Faibussowitsch     PetscCall(SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs));
462ab8d36c9SPeter Brune   } else {
4639566063dSJacob Faibussowitsch     PetscCall(PetscSNPrintf(tprefix,sizeof(tprefix),"fas_levels_%d_",(int)fas->level));
4649566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, optionsprefix));
4659566063dSJacob Faibussowitsch     PetscCall(SNESAppendOptionsPrefix(nsmooth, tprefix));
4669566063dSJacob Faibussowitsch     PetscCall(SNESSetType(nsmooth, SNESNRICHARDSON));
4679566063dSJacob Faibussowitsch     PetscCall(SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs));
468ab8d36c9SPeter Brune   }
4699566063dSJacob Faibussowitsch   PetscCall(PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1));
4709566063dSJacob Faibussowitsch   PetscCall(PetscLogObjectParent((PetscObject)snes,(PetscObject)nsmooth));
4719566063dSJacob Faibussowitsch   PetscCall(PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth));
4729566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposedDataSetInt((PetscObject) nsmooth, PetscMGLevelId, fas->level));
473ab8d36c9SPeter Brune   *smooth = nsmooth;
474ab8d36c9SPeter Brune   PetscFunctionReturn(0);
475ab8d36c9SPeter Brune }
476ab8d36c9SPeter Brune 
477ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */
478ab8d36c9SPeter Brune 
479ab8d36c9SPeter Brune /*@
480ab8d36c9SPeter Brune    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
481ab8d36c9SPeter Brune 
482ab8d36c9SPeter Brune    Logically Collective on SNES
483ab8d36c9SPeter Brune 
484ab8d36c9SPeter Brune    Input Parameters:
485ab8d36c9SPeter Brune +  snes   - the multigrid context
486ab8d36c9SPeter Brune .  level  - the level to set the number of cycles on
487ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
488ab8d36c9SPeter Brune 
489ab8d36c9SPeter Brune    Level: advanced
490ab8d36c9SPeter Brune 
491ab8d36c9SPeter Brune .seealso: SNESFASSetCycles()
492ab8d36c9SPeter Brune @*/
49322d28d08SBarry Smith PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
49422d28d08SBarry Smith {
495f833ba53SLisandro Dalcin   SNES_FAS       *fas;
49622d28d08SBarry Smith 
497ab8d36c9SPeter Brune   PetscFunctionBegin;
498f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
499f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
500ab8d36c9SPeter Brune   fas->n_cycles = cycles;
5019566063dSJacob Faibussowitsch   PetscCall(SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs));
502ab8d36c9SPeter Brune   PetscFunctionReturn(0);
503ab8d36c9SPeter Brune }
504ab8d36c9SPeter Brune 
505ab8d36c9SPeter Brune /*@
506ab8d36c9SPeter Brune    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
507ab8d36c9SPeter Brune 
508ab8d36c9SPeter Brune    Logically Collective on SNES
509ab8d36c9SPeter Brune 
510ab8d36c9SPeter Brune    Input Parameters:
511ab8d36c9SPeter Brune .  snes   - the multigrid context
512ab8d36c9SPeter Brune 
513ab8d36c9SPeter Brune    Output Parameters:
514ab8d36c9SPeter Brune .  smooth - the smoother
515ab8d36c9SPeter Brune 
516ab8d36c9SPeter Brune    Level: advanced
517ab8d36c9SPeter Brune 
518ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
519ab8d36c9SPeter Brune @*/
520ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
521ab8d36c9SPeter Brune {
522ab8d36c9SPeter Brune   SNES_FAS *fas;
5235fd66863SKarl Rupp 
524ab8d36c9SPeter Brune   PetscFunctionBegin;
525f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
526f833ba53SLisandro Dalcin   PetscValidPointer(smooth,2);
527ab8d36c9SPeter Brune   fas     = (SNES_FAS*)snes->data;
528ab8d36c9SPeter Brune   *smooth = fas->smoothd;
529ab8d36c9SPeter Brune   PetscFunctionReturn(0);
530ab8d36c9SPeter Brune }
531ab8d36c9SPeter Brune /*@
532ab8d36c9SPeter Brune    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
533ab8d36c9SPeter Brune 
534ab8d36c9SPeter Brune    Logically Collective on SNES
535ab8d36c9SPeter Brune 
536ab8d36c9SPeter Brune    Input Parameters:
537ab8d36c9SPeter Brune .  snes   - the multigrid context
538ab8d36c9SPeter Brune 
539ab8d36c9SPeter Brune    Output Parameters:
540ab8d36c9SPeter Brune .  smoothu - the smoother
541ab8d36c9SPeter Brune 
542ab8d36c9SPeter Brune    Notes:
543ab8d36c9SPeter Brune    Returns the downsmoother if no up smoother is available.  This enables transparent
544ab8d36c9SPeter Brune    default behavior in the process of the solve.
545ab8d36c9SPeter Brune 
546ab8d36c9SPeter Brune    Level: advanced
547ab8d36c9SPeter Brune 
548ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
549ab8d36c9SPeter Brune @*/
550ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
551ab8d36c9SPeter Brune {
552ab8d36c9SPeter Brune   SNES_FAS *fas;
5535fd66863SKarl Rupp 
554ab8d36c9SPeter Brune   PetscFunctionBegin;
555f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
556f833ba53SLisandro Dalcin   PetscValidPointer(smoothu,2);
557ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
5581aa26658SKarl Rupp   if (!fas->smoothu) *smoothu = fas->smoothd;
5591aa26658SKarl Rupp   else *smoothu = fas->smoothu;
560ab8d36c9SPeter Brune   PetscFunctionReturn(0);
561ab8d36c9SPeter Brune }
562ab8d36c9SPeter Brune 
563ab8d36c9SPeter Brune /*@
564ab8d36c9SPeter Brune    SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
565ab8d36c9SPeter Brune 
566ab8d36c9SPeter Brune    Logically Collective on SNES
567ab8d36c9SPeter Brune 
568ab8d36c9SPeter Brune    Input Parameters:
569ab8d36c9SPeter Brune .  snes   - the multigrid context
570ab8d36c9SPeter Brune 
571ab8d36c9SPeter Brune    Output Parameters:
572ab8d36c9SPeter Brune .  smoothd - the smoother
573ab8d36c9SPeter Brune 
574ab8d36c9SPeter Brune    Level: advanced
575ab8d36c9SPeter Brune 
576ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
577ab8d36c9SPeter Brune @*/
578ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
579ab8d36c9SPeter Brune {
580ab8d36c9SPeter Brune   SNES_FAS *fas;
5815fd66863SKarl Rupp 
582ab8d36c9SPeter Brune   PetscFunctionBegin;
583f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
584f833ba53SLisandro Dalcin   PetscValidPointer(smoothd,2);
585ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
586ab8d36c9SPeter Brune   *smoothd = fas->smoothd;
587ab8d36c9SPeter Brune   PetscFunctionReturn(0);
588ab8d36c9SPeter Brune }
589ab8d36c9SPeter Brune 
590ab8d36c9SPeter Brune /*@
591ab8d36c9SPeter Brune    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level
592ab8d36c9SPeter Brune 
593ab8d36c9SPeter Brune    Logically Collective on SNES
594ab8d36c9SPeter Brune 
595ab8d36c9SPeter Brune    Input Parameters:
596ab8d36c9SPeter Brune .  snes   - the multigrid context
597ab8d36c9SPeter Brune 
598ab8d36c9SPeter Brune    Output Parameters:
599ab8d36c9SPeter Brune .  correction - the coarse correction on this level
600ab8d36c9SPeter Brune 
601ab8d36c9SPeter Brune    Notes:
6020298fd71SBarry Smith    Returns NULL on the coarsest level.
603ab8d36c9SPeter Brune 
604ab8d36c9SPeter Brune    Level: advanced
605ab8d36c9SPeter Brune 
606ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
607ab8d36c9SPeter Brune @*/
608ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
609ab8d36c9SPeter Brune {
610ab8d36c9SPeter Brune   SNES_FAS *fas;
6115fd66863SKarl Rupp 
612ab8d36c9SPeter Brune   PetscFunctionBegin;
613f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
614f833ba53SLisandro Dalcin   PetscValidPointer(correction,2);
615ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
616ab8d36c9SPeter Brune   *correction = fas->next;
617ab8d36c9SPeter Brune   PetscFunctionReturn(0);
618ab8d36c9SPeter Brune }
619ab8d36c9SPeter Brune 
620ab8d36c9SPeter Brune /*@
621ab8d36c9SPeter Brune    SNESFASCycleGetInterpolation - Gets the interpolation on this level
622ab8d36c9SPeter Brune 
623ab8d36c9SPeter Brune    Logically Collective on SNES
624ab8d36c9SPeter Brune 
625ab8d36c9SPeter Brune    Input Parameters:
626ab8d36c9SPeter Brune .  snes   - the multigrid context
627ab8d36c9SPeter Brune 
628ab8d36c9SPeter Brune    Output Parameters:
629ab8d36c9SPeter Brune .  mat    - the interpolation operator on this level
630ab8d36c9SPeter Brune 
631ab8d36c9SPeter Brune    Level: developer
632ab8d36c9SPeter Brune 
633ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
634ab8d36c9SPeter Brune @*/
635ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
636ab8d36c9SPeter Brune {
637ab8d36c9SPeter Brune   SNES_FAS *fas;
6385fd66863SKarl Rupp 
639ab8d36c9SPeter Brune   PetscFunctionBegin;
640f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
641f833ba53SLisandro Dalcin   PetscValidPointer(mat,2);
642ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
643ab8d36c9SPeter Brune   *mat = fas->interpolate;
644ab8d36c9SPeter Brune   PetscFunctionReturn(0);
645ab8d36c9SPeter Brune }
646ab8d36c9SPeter Brune 
647ab8d36c9SPeter Brune /*@
648ab8d36c9SPeter Brune    SNESFASCycleGetRestriction - Gets the restriction on this level
649ab8d36c9SPeter Brune 
650ab8d36c9SPeter Brune    Logically Collective on SNES
651ab8d36c9SPeter Brune 
652ab8d36c9SPeter Brune    Input Parameters:
653ab8d36c9SPeter Brune .  snes   - the multigrid context
654ab8d36c9SPeter Brune 
655ab8d36c9SPeter Brune    Output Parameters:
656ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
657ab8d36c9SPeter Brune 
658ab8d36c9SPeter Brune    Level: developer
659ab8d36c9SPeter Brune 
660ab8d36c9SPeter Brune .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
661ab8d36c9SPeter Brune @*/
662ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
663ab8d36c9SPeter Brune {
664ab8d36c9SPeter Brune   SNES_FAS *fas;
6655fd66863SKarl Rupp 
666ab8d36c9SPeter Brune   PetscFunctionBegin;
667f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
668f833ba53SLisandro Dalcin   PetscValidPointer(mat,2);
669ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
670ab8d36c9SPeter Brune   *mat = fas->restrct;
671ab8d36c9SPeter Brune   PetscFunctionReturn(0);
672ab8d36c9SPeter Brune }
673ab8d36c9SPeter Brune 
674ab8d36c9SPeter Brune /*@
675ab8d36c9SPeter Brune    SNESFASCycleGetInjection - Gets the injection on this level
676ab8d36c9SPeter Brune 
677ab8d36c9SPeter Brune    Logically Collective on SNES
678ab8d36c9SPeter Brune 
679ab8d36c9SPeter Brune    Input Parameters:
680ab8d36c9SPeter Brune .  snes   - the multigrid context
681ab8d36c9SPeter Brune 
682ab8d36c9SPeter Brune    Output Parameters:
683ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
684ab8d36c9SPeter Brune 
685ab8d36c9SPeter Brune    Level: developer
686ab8d36c9SPeter Brune 
687ab8d36c9SPeter Brune .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
688ab8d36c9SPeter Brune @*/
689ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
690ab8d36c9SPeter Brune {
691ab8d36c9SPeter Brune   SNES_FAS *fas;
6925fd66863SKarl Rupp 
693ab8d36c9SPeter Brune   PetscFunctionBegin;
694f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
695f833ba53SLisandro Dalcin   PetscValidPointer(mat,2);
696ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
697ab8d36c9SPeter Brune   *mat = fas->inject;
698ab8d36c9SPeter Brune   PetscFunctionReturn(0);
699ab8d36c9SPeter Brune }
700ab8d36c9SPeter Brune 
701ab8d36c9SPeter Brune /*@
702ab8d36c9SPeter Brune    SNESFASCycleGetRScale - Gets the injection on this level
703ab8d36c9SPeter Brune 
704ab8d36c9SPeter Brune    Logically Collective on SNES
705ab8d36c9SPeter Brune 
706ab8d36c9SPeter Brune    Input Parameters:
707ab8d36c9SPeter Brune .  snes   - the multigrid context
708ab8d36c9SPeter Brune 
709ab8d36c9SPeter Brune    Output Parameters:
710ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
711ab8d36c9SPeter Brune 
712ab8d36c9SPeter Brune    Level: developer
713ab8d36c9SPeter Brune 
714ab8d36c9SPeter Brune .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
715ab8d36c9SPeter Brune @*/
716ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
717ab8d36c9SPeter Brune {
718ab8d36c9SPeter Brune   SNES_FAS *fas;
7195fd66863SKarl Rupp 
720ab8d36c9SPeter Brune   PetscFunctionBegin;
721f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
722f833ba53SLisandro Dalcin   PetscValidPointer(vec,2);
723ab8d36c9SPeter Brune   fas  = (SNES_FAS*)snes->data;
724ab8d36c9SPeter Brune   *vec = fas->rscale;
725ab8d36c9SPeter Brune   PetscFunctionReturn(0);
726ab8d36c9SPeter Brune }
727ab8d36c9SPeter Brune 
728ab8d36c9SPeter Brune /*@
729ab8d36c9SPeter Brune    SNESFASCycleIsFine - Determines if a given cycle is the fine level.
730ab8d36c9SPeter Brune 
731ab8d36c9SPeter Brune    Logically Collective on SNES
732ab8d36c9SPeter Brune 
733ab8d36c9SPeter Brune    Input Parameters:
734ab8d36c9SPeter Brune .  snes   - the FAS context
735ab8d36c9SPeter Brune 
736ab8d36c9SPeter Brune    Output Parameters:
737ab8d36c9SPeter Brune .  flg - indicates if this is the fine level or not
738ab8d36c9SPeter Brune 
739ab8d36c9SPeter Brune    Level: advanced
740ab8d36c9SPeter Brune 
741ab8d36c9SPeter Brune .seealso: SNESFASSetLevels()
742ab8d36c9SPeter Brune @*/
743ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
744ab8d36c9SPeter Brune {
745ab8d36c9SPeter Brune   SNES_FAS *fas;
7465fd66863SKarl Rupp 
747ab8d36c9SPeter Brune   PetscFunctionBegin;
748f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
749534a8f05SLisandro Dalcin   PetscValidBoolPointer(flg,2);
750ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
7511aa26658SKarl Rupp   if (fas->level == fas->levels - 1) *flg = PETSC_TRUE;
7521aa26658SKarl Rupp   else *flg = PETSC_FALSE;
753ab8d36c9SPeter Brune   PetscFunctionReturn(0);
754ab8d36c9SPeter Brune }
755ab8d36c9SPeter Brune 
756ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */
757ab8d36c9SPeter Brune 
758ab8d36c9SPeter Brune /*@
759ab8d36c9SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
760ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
761ab8d36c9SPeter Brune 
762ab8d36c9SPeter Brune    Input Parameters:
763ab8d36c9SPeter Brune +  snes      - the multigrid context
764ab8d36c9SPeter Brune .  mat       - the interpolation operator
765ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
766ab8d36c9SPeter Brune 
767ab8d36c9SPeter Brune    Level: advanced
768ab8d36c9SPeter Brune 
769ab8d36c9SPeter Brune    Notes:
770ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the restriction
771ab8d36c9SPeter Brune     for the same level.
772ab8d36c9SPeter Brune 
773ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
774ab8d36c9SPeter Brune     out from the matrix size which one it is.
775ab8d36c9SPeter Brune 
776ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
777ab8d36c9SPeter Brune @*/
7780adebc6cSBarry Smith PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
7790adebc6cSBarry Smith {
78022d28d08SBarry Smith   SNES_FAS       *fas;
781ab8d36c9SPeter Brune   SNES           levelsnes;
78222d28d08SBarry Smith 
783ab8d36c9SPeter Brune   PetscFunctionBegin;
784f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
785f833ba53SLisandro Dalcin   if (mat) PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
7869566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
787ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
7889566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
7899566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->interpolate));
790ab8d36c9SPeter Brune   fas->interpolate = mat;
791ab8d36c9SPeter Brune   PetscFunctionReturn(0);
792ab8d36c9SPeter Brune }
793ab8d36c9SPeter Brune 
794ab8d36c9SPeter Brune /*@
795ab8d36c9SPeter Brune    SNESFASGetInterpolation - Gets the matrix used to calculate the
796ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
797ab8d36c9SPeter Brune 
798ab8d36c9SPeter Brune    Input Parameters:
799ab8d36c9SPeter Brune +  snes      - the multigrid context
800ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
801ab8d36c9SPeter Brune 
802ab8d36c9SPeter Brune    Output Parameters:
803ab8d36c9SPeter Brune .  mat       - the interpolation operator
804ab8d36c9SPeter Brune 
805ab8d36c9SPeter Brune    Level: advanced
806ab8d36c9SPeter Brune 
807ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
808ab8d36c9SPeter Brune @*/
80922d28d08SBarry Smith PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
81022d28d08SBarry Smith {
81122d28d08SBarry Smith   SNES_FAS       *fas;
812ab8d36c9SPeter Brune   SNES           levelsnes;
81322d28d08SBarry Smith 
814ab8d36c9SPeter Brune   PetscFunctionBegin;
815f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
816f833ba53SLisandro Dalcin   PetscValidPointer(mat,3);
8179566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
818ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
819ab8d36c9SPeter Brune   *mat = fas->interpolate;
820ab8d36c9SPeter Brune   PetscFunctionReturn(0);
821ab8d36c9SPeter Brune }
822ab8d36c9SPeter Brune 
823ab8d36c9SPeter Brune /*@
824ab8d36c9SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
825ab8d36c9SPeter Brune    from level l to l-1.
826ab8d36c9SPeter Brune 
827ab8d36c9SPeter Brune    Input Parameters:
828ab8d36c9SPeter Brune +  snes  - the multigrid context
829ab8d36c9SPeter Brune .  mat   - the restriction matrix
830ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
831ab8d36c9SPeter Brune 
832ab8d36c9SPeter Brune    Level: advanced
833ab8d36c9SPeter Brune 
834ab8d36c9SPeter Brune    Notes:
835ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the interpolation
836ab8d36c9SPeter Brune     for the same level.
837ab8d36c9SPeter Brune 
838ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
839ab8d36c9SPeter Brune     out from the matrix size which one it is.
840ab8d36c9SPeter Brune 
841ab8d36c9SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
842ab8d36c9SPeter Brune     is used.
843ab8d36c9SPeter Brune 
844ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
845ab8d36c9SPeter Brune @*/
84622d28d08SBarry Smith PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
84722d28d08SBarry Smith {
84822d28d08SBarry Smith   SNES_FAS       *fas;
849ab8d36c9SPeter Brune   SNES           levelsnes;
85022d28d08SBarry Smith 
851ab8d36c9SPeter Brune   PetscFunctionBegin;
852f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
853f833ba53SLisandro Dalcin   if (mat) PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
8549566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
855ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
8569566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
8579566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->restrct));
858ab8d36c9SPeter Brune   fas->restrct = mat;
859ab8d36c9SPeter Brune   PetscFunctionReturn(0);
860ab8d36c9SPeter Brune }
861ab8d36c9SPeter Brune 
862ab8d36c9SPeter Brune /*@
863ab8d36c9SPeter Brune    SNESFASGetRestriction - Gets the matrix used to calculate the
864ab8d36c9SPeter Brune    restriction from l to the l-1th level
865ab8d36c9SPeter Brune 
866ab8d36c9SPeter Brune    Input Parameters:
867ab8d36c9SPeter Brune +  snes      - the multigrid context
868ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
869ab8d36c9SPeter Brune 
870ab8d36c9SPeter Brune    Output Parameters:
871ab8d36c9SPeter Brune .  mat       - the interpolation operator
872ab8d36c9SPeter Brune 
873ab8d36c9SPeter Brune    Level: advanced
874ab8d36c9SPeter Brune 
875ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
876ab8d36c9SPeter Brune @*/
87722d28d08SBarry Smith PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
87822d28d08SBarry Smith {
87922d28d08SBarry Smith   SNES_FAS       *fas;
880ab8d36c9SPeter Brune   SNES           levelsnes;
88122d28d08SBarry Smith 
882ab8d36c9SPeter Brune   PetscFunctionBegin;
883f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
884f833ba53SLisandro Dalcin   PetscValidPointer(mat,3);
8859566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
886ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
887ab8d36c9SPeter Brune   *mat = fas->restrct;
888ab8d36c9SPeter Brune   PetscFunctionReturn(0);
889ab8d36c9SPeter Brune }
890ab8d36c9SPeter Brune 
891ab8d36c9SPeter Brune /*@
892ab8d36c9SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
893ab8d36c9SPeter Brune    from level l to l-1.
894ab8d36c9SPeter Brune 
895ab8d36c9SPeter Brune    Input Parameters:
896ab8d36c9SPeter Brune  +  snes  - the multigrid context
897ab8d36c9SPeter Brune .  mat   - the restriction matrix
898ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
899ab8d36c9SPeter Brune 
900ab8d36c9SPeter Brune    Level: advanced
901ab8d36c9SPeter Brune 
902ab8d36c9SPeter Brune    Notes:
903ab8d36c9SPeter Brune          If you do not set this, the restriction and rscale is used to
904ab8d36c9SPeter Brune    project the solution instead.
905ab8d36c9SPeter Brune 
906ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
907ab8d36c9SPeter Brune @*/
90822d28d08SBarry Smith PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
90922d28d08SBarry Smith {
91022d28d08SBarry Smith   SNES_FAS       *fas;
911ab8d36c9SPeter Brune   SNES           levelsnes;
91222d28d08SBarry Smith 
913ab8d36c9SPeter Brune   PetscFunctionBegin;
914f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
915f833ba53SLisandro Dalcin   if (mat) PetscValidHeaderSpecific(mat,MAT_CLASSID,3);
9169566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
917ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
9189566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)mat));
9199566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&fas->inject));
9201aa26658SKarl Rupp 
921ab8d36c9SPeter Brune   fas->inject = mat;
922ab8d36c9SPeter Brune   PetscFunctionReturn(0);
923ab8d36c9SPeter Brune }
924ab8d36c9SPeter Brune 
925ab8d36c9SPeter Brune /*@
926ab8d36c9SPeter Brune    SNESFASGetInjection - Gets the matrix used to calculate the
927ab8d36c9SPeter Brune    injection from l-1 to the lth level
928ab8d36c9SPeter Brune 
929ab8d36c9SPeter Brune    Input Parameters:
930ab8d36c9SPeter Brune +  snes      - the multigrid context
931ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
932ab8d36c9SPeter Brune 
933ab8d36c9SPeter Brune    Output Parameters:
934ab8d36c9SPeter Brune .  mat       - the injection operator
935ab8d36c9SPeter Brune 
936ab8d36c9SPeter Brune    Level: advanced
937ab8d36c9SPeter Brune 
938ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
939ab8d36c9SPeter Brune @*/
94022d28d08SBarry Smith PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
94122d28d08SBarry Smith {
94222d28d08SBarry Smith   SNES_FAS       *fas;
943ab8d36c9SPeter Brune   SNES           levelsnes;
94422d28d08SBarry Smith 
945ab8d36c9SPeter Brune   PetscFunctionBegin;
946f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
947f833ba53SLisandro Dalcin   PetscValidPointer(mat,3);
9489566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
949ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
950ab8d36c9SPeter Brune   *mat = fas->inject;
951ab8d36c9SPeter Brune   PetscFunctionReturn(0);
952ab8d36c9SPeter Brune }
953ab8d36c9SPeter Brune 
954ab8d36c9SPeter Brune /*@
955ab8d36c9SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
956ab8d36c9SPeter Brune    operator from level l to l-1.
957ab8d36c9SPeter Brune 
958ab8d36c9SPeter Brune    Input Parameters:
959ab8d36c9SPeter Brune +  snes   - the multigrid context
960ab8d36c9SPeter Brune .  rscale - the restriction scaling
961ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
962ab8d36c9SPeter Brune 
963ab8d36c9SPeter Brune    Level: advanced
964ab8d36c9SPeter Brune 
965ab8d36c9SPeter Brune    Notes:
966ab8d36c9SPeter Brune          This is only used in the case that the injection is not set.
967ab8d36c9SPeter Brune 
968ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
969ab8d36c9SPeter Brune @*/
97022d28d08SBarry Smith PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
97122d28d08SBarry Smith {
972ab8d36c9SPeter Brune   SNES_FAS       *fas;
973ab8d36c9SPeter Brune   SNES           levelsnes;
97422d28d08SBarry Smith 
975ab8d36c9SPeter Brune   PetscFunctionBegin;
976f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
977f833ba53SLisandro Dalcin   if (rscale) PetscValidHeaderSpecific(rscale,VEC_CLASSID,3);
9789566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
979ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
9809566063dSJacob Faibussowitsch   PetscCall(PetscObjectReference((PetscObject)rscale));
9819566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&fas->rscale));
982ab8d36c9SPeter Brune   fas->rscale = rscale;
983ab8d36c9SPeter Brune   PetscFunctionReturn(0);
984ab8d36c9SPeter Brune }
985ab8d36c9SPeter Brune 
986ab8d36c9SPeter Brune /*@
987ab8d36c9SPeter Brune    SNESFASGetSmoother - Gets the default smoother on a level.
988ab8d36c9SPeter Brune 
989ab8d36c9SPeter Brune    Input Parameters:
990ab8d36c9SPeter Brune +  snes   - the multigrid context
991ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
992ab8d36c9SPeter Brune 
993ab8d36c9SPeter Brune    Output Parameters:
994ab8d36c9SPeter Brune    smooth  - the smoother
995ab8d36c9SPeter Brune 
996ab8d36c9SPeter Brune    Level: advanced
997ab8d36c9SPeter Brune 
998ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
999ab8d36c9SPeter Brune @*/
100022d28d08SBarry Smith PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
100122d28d08SBarry Smith {
1002ab8d36c9SPeter Brune   SNES_FAS       *fas;
1003ab8d36c9SPeter Brune   SNES           levelsnes;
100422d28d08SBarry Smith 
1005ab8d36c9SPeter Brune   PetscFunctionBegin;
1006f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
1007f833ba53SLisandro Dalcin   PetscValidPointer(smooth,3);
10089566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1009ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1010ab8d36c9SPeter Brune   if (!fas->smoothd) {
10119566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1012ab8d36c9SPeter Brune   }
1013ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1014ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1015ab8d36c9SPeter Brune }
1016ab8d36c9SPeter Brune 
1017ab8d36c9SPeter Brune /*@
1018ab8d36c9SPeter Brune    SNESFASGetSmootherDown - Gets the downsmoother on a level.
1019ab8d36c9SPeter Brune 
1020ab8d36c9SPeter Brune    Input Parameters:
1021ab8d36c9SPeter Brune +  snes   - the multigrid context
1022ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
1023ab8d36c9SPeter Brune 
1024ab8d36c9SPeter Brune    Output Parameters:
1025ab8d36c9SPeter Brune    smooth  - the smoother
1026ab8d36c9SPeter Brune 
1027ab8d36c9SPeter Brune    Level: advanced
1028ab8d36c9SPeter Brune 
1029ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1030ab8d36c9SPeter Brune @*/
103122d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
103222d28d08SBarry Smith {
1033ab8d36c9SPeter Brune   SNES_FAS       *fas;
1034ab8d36c9SPeter Brune   SNES           levelsnes;
103522d28d08SBarry Smith 
1036ab8d36c9SPeter Brune   PetscFunctionBegin;
1037f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
1038f833ba53SLisandro Dalcin   PetscValidPointer(smooth,3);
10399566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1040ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1041ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1042ab8d36c9SPeter Brune   if (!fas->smoothd) {
10439566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1044ab8d36c9SPeter Brune   }
1045ab8d36c9SPeter Brune   if (!fas->smoothu) {
10469566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
1047ab8d36c9SPeter Brune   }
1048ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1049ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1050ab8d36c9SPeter Brune }
1051ab8d36c9SPeter Brune 
1052ab8d36c9SPeter Brune /*@
1053ab8d36c9SPeter Brune    SNESFASGetSmootherUp - Gets the upsmoother on a level.
1054ab8d36c9SPeter Brune 
1055ab8d36c9SPeter Brune    Input Parameters:
1056ab8d36c9SPeter Brune +  snes   - the multigrid context
1057ab8d36c9SPeter Brune -  level  - the level (0 is coarsest)
1058ab8d36c9SPeter Brune 
1059ab8d36c9SPeter Brune    Output Parameters:
1060ab8d36c9SPeter Brune    smooth  - the smoother
1061ab8d36c9SPeter Brune 
1062ab8d36c9SPeter Brune    Level: advanced
1063ab8d36c9SPeter Brune 
1064ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1065ab8d36c9SPeter Brune @*/
106622d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
106722d28d08SBarry Smith {
1068ab8d36c9SPeter Brune   SNES_FAS       *fas;
1069ab8d36c9SPeter Brune   SNES           levelsnes;
107022d28d08SBarry Smith 
1071ab8d36c9SPeter Brune   PetscFunctionBegin;
1072f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
1073f833ba53SLisandro Dalcin   PetscValidPointer(smooth,3);
10749566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, level, &levelsnes));
1075ab8d36c9SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1076ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1077ab8d36c9SPeter Brune   if (!fas->smoothd) {
10789566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1079ab8d36c9SPeter Brune   }
1080ab8d36c9SPeter Brune   if (!fas->smoothu) {
10819566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothu));
1082ab8d36c9SPeter Brune   }
1083ab8d36c9SPeter Brune   *smooth = fas->smoothu;
1084ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1085ab8d36c9SPeter Brune }
1086d6ad1212SPeter Brune 
1087d6ad1212SPeter Brune /*@
1088d6ad1212SPeter Brune   SNESFASGetCoarseSolve - Gets the coarsest solver.
1089d6ad1212SPeter Brune 
1090d6ad1212SPeter Brune   Input Parameters:
1091a3a80b83SMatthew G. Knepley . snes - the multigrid context
1092d6ad1212SPeter Brune 
1093d6ad1212SPeter Brune   Output Parameters:
1094a3a80b83SMatthew G. Knepley . coarse - the coarse-level solver
1095d6ad1212SPeter Brune 
1096d6ad1212SPeter Brune   Level: advanced
1097d6ad1212SPeter Brune 
1098d6ad1212SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1099d6ad1212SPeter Brune @*/
1100a3a80b83SMatthew G. Knepley PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *coarse)
110122d28d08SBarry Smith {
1102d6ad1212SPeter Brune   SNES_FAS       *fas;
1103d6ad1212SPeter Brune   SNES           levelsnes;
110422d28d08SBarry Smith 
1105d6ad1212SPeter Brune   PetscFunctionBegin;
1106f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
1107064a246eSJacob Faibussowitsch   PetscValidPointer(coarse,2);
11089566063dSJacob Faibussowitsch   PetscCall(SNESFASGetCycleSNES(snes, 0, &levelsnes));
1109d6ad1212SPeter Brune   fas  = (SNES_FAS*)levelsnes->data;
1110d6ad1212SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1111d6ad1212SPeter Brune   if (!fas->smoothd) {
11129566063dSJacob Faibussowitsch     PetscCall(SNESFASCycleCreateSmoother_Private(levelsnes, &fas->smoothd));
1113d6ad1212SPeter Brune   }
1114a3a80b83SMatthew G. Knepley   *coarse = fas->smoothd;
1115d6ad1212SPeter Brune   PetscFunctionReturn(0);
1116d6ad1212SPeter Brune }
1117928e959bSPeter Brune 
1118928e959bSPeter Brune /*@
1119928e959bSPeter Brune    SNESFASFullSetDownSweep - Smooth during the initial downsweep for SNESFAS
1120928e959bSPeter Brune 
1121928e959bSPeter Brune    Logically Collective on SNES
1122928e959bSPeter Brune 
1123928e959bSPeter Brune    Input Parameters:
1124928e959bSPeter Brune +  snes - the multigrid context
1125928e959bSPeter Brune -  swp - whether to downsweep or not
1126928e959bSPeter Brune 
1127928e959bSPeter Brune    Options Database Key:
1128928e959bSPeter Brune .  -snes_fas_full_downsweep - Sets number of pre-smoothing steps
1129928e959bSPeter Brune 
1130928e959bSPeter Brune    Level: advanced
1131928e959bSPeter Brune 
1132928e959bSPeter Brune .seealso: SNESFASSetNumberSmoothUp()
1133928e959bSPeter Brune @*/
1134928e959bSPeter Brune PetscErrorCode SNESFASFullSetDownSweep(SNES snes,PetscBool swp)
1135928e959bSPeter Brune {
1136f833ba53SLisandro Dalcin   SNES_FAS       *fas;
1137928e959bSPeter Brune 
1138928e959bSPeter Brune   PetscFunctionBegin;
1139f833ba53SLisandro Dalcin   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
1140f833ba53SLisandro Dalcin   fas = (SNES_FAS*)snes->data;
1141928e959bSPeter Brune   fas->full_downsweep = swp;
1142928e959bSPeter Brune   if (fas->next) {
11439566063dSJacob Faibussowitsch     PetscCall(SNESFASFullSetDownSweep(fas->next,swp));
1144928e959bSPeter Brune   }
1145928e959bSPeter Brune   PetscFunctionReturn(0);
1146928e959bSPeter Brune }
11476dfbafefSToby Isaac 
11486dfbafefSToby Isaac /*@
11496dfbafefSToby Isaac    SNESFASFullSetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles
11506dfbafefSToby Isaac 
11516dfbafefSToby Isaac    Logically Collective on SNES
11526dfbafefSToby Isaac 
11536dfbafefSToby Isaac    Input Parameters:
11546dfbafefSToby Isaac +  snes - the multigrid context
11556dfbafefSToby Isaac -  total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
11566dfbafefSToby Isaac 
11576dfbafefSToby Isaac    Options Database Key:
115867b8a455SSatish Balay .  -snes_fas_full_total - Use total restriction and interpolation on the initial down and up sweeps for the full FAS cycle
11596dfbafefSToby Isaac 
11606dfbafefSToby Isaac    Level: advanced
11616dfbafefSToby Isaac 
11626dfbafefSToby Isaac    Note: this option is only significant if the interpolation of a coarse correction (MatInterpolate()) is significantly different from total solution interpolation (DMInterpolateSolution()).
11636dfbafefSToby Isaac 
11646dfbafefSToby Isaac .seealso: SNESFASSetNumberSmoothUp(), DMInterpolateSolution()
11656dfbafefSToby Isaac @*/
11666dfbafefSToby Isaac PetscErrorCode SNESFASFullSetTotal(SNES snes,PetscBool total)
11676dfbafefSToby Isaac {
11686dfbafefSToby Isaac   SNES_FAS       *fas;
11696dfbafefSToby Isaac 
11706dfbafefSToby Isaac   PetscFunctionBegin;
11716dfbafefSToby Isaac   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
11726dfbafefSToby Isaac   fas = (SNES_FAS*)snes->data;
11736dfbafefSToby Isaac   fas->full_total = total;
11746dfbafefSToby Isaac   if (fas->next) {
11759566063dSJacob Faibussowitsch     PetscCall(SNESFASFullSetTotal(fas->next,total));
11766dfbafefSToby Isaac   }
11776dfbafefSToby Isaac   PetscFunctionReturn(0);
11786dfbafefSToby Isaac }
11796dfbafefSToby Isaac 
11806dfbafefSToby Isaac /*@
11816dfbafefSToby Isaac    SNESFASFullGetTotal - Use total residual restriction and total interpolation on the initial down and up sweep of full FAS cycles
11826dfbafefSToby Isaac 
11836dfbafefSToby Isaac    Logically Collective on SNES
11846dfbafefSToby Isaac 
11856dfbafefSToby Isaac    Input Parameters:
11866dfbafefSToby Isaac .  snes - the multigrid context
11876dfbafefSToby Isaac 
11886dfbafefSToby Isaac    Output:
11896dfbafefSToby Isaac .  total - whether to use total restriction / interpolatiaon or not (the alternative is defect restriction and correction interpolation)
11906dfbafefSToby Isaac 
11916dfbafefSToby Isaac    Level: advanced
11926dfbafefSToby Isaac 
11936dfbafefSToby Isaac .seealso: SNESFASSetNumberSmoothUp(), DMInterpolateSolution(), SNESFullSetTotal()
11946dfbafefSToby Isaac @*/
11956dfbafefSToby Isaac PetscErrorCode SNESFASFullGetTotal(SNES snes,PetscBool *total)
11966dfbafefSToby Isaac {
11976dfbafefSToby Isaac   SNES_FAS       *fas;
11986dfbafefSToby Isaac 
11996dfbafefSToby Isaac   PetscFunctionBegin;
12006dfbafefSToby Isaac   PetscValidHeaderSpecificType(snes,SNES_CLASSID,1,SNESFAS);
12016dfbafefSToby Isaac   fas = (SNES_FAS*)snes->data;
12026dfbafefSToby Isaac   *total = fas->full_total;
12036dfbafefSToby Isaac   PetscFunctionReturn(0);
12046dfbafefSToby Isaac }
1205