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