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