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