xref: /petsc/src/snes/impls/fas/fasfunc.c (revision 5fd668637986a8d8518383a9159eebc368e1d5b4)
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;
158*5fd66863SKarl Rupp 
159ab8d36c9SPeter Brune   PetscFunctionBegin;
160ab8d36c9SPeter Brune   *levels = fas->levels;
161ab8d36c9SPeter Brune   PetscFunctionReturn(0);
162ab8d36c9SPeter Brune }
163ab8d36c9SPeter Brune 
164ab8d36c9SPeter Brune 
165ab8d36c9SPeter Brune #undef __FUNCT__
166ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetCycleSNES"
167ab8d36c9SPeter Brune /*@
168ab8d36c9SPeter Brune    SNESFASGetCycleSNES - Gets the SNES corresponding to a particular
169ab8d36c9SPeter Brune    level of the FAS hierarchy.
170ab8d36c9SPeter Brune 
171ab8d36c9SPeter Brune    Input Parameters:
172ab8d36c9SPeter Brune +  snes    - the multigrid context
173ab8d36c9SPeter Brune    level   - the level to get
174ab8d36c9SPeter Brune -  lsnes   - whether to use the nonlinear smoother or not
175ab8d36c9SPeter Brune 
176ab8d36c9SPeter Brune    Level: advanced
177ab8d36c9SPeter Brune 
178ab8d36c9SPeter Brune .keywords: FAS, MG, set, cycles, Gauss-Seidel, multigrid
179ab8d36c9SPeter Brune 
180ab8d36c9SPeter Brune .seealso: SNESFASSetLevels(), SNESFASGetLevels()
181ab8d36c9SPeter Brune @*/
18222d28d08SBarry Smith PetscErrorCode SNESFASGetCycleSNES(SNES snes,PetscInt level,SNES *lsnes)
18322d28d08SBarry Smith {
184ab8d36c9SPeter Brune   SNES_FAS *fas = (SNES_FAS*)snes->data;
185ab8d36c9SPeter Brune   PetscInt i;
186ab8d36c9SPeter Brune 
187ab8d36c9SPeter Brune   PetscFunctionBegin;
188ab8d36c9SPeter 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);
18922d28d08SBarry 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);
190ab8d36c9SPeter Brune 
191ab8d36c9SPeter Brune   *lsnes = snes;
192ab8d36c9SPeter Brune   for (i = fas->level; i > level; i--) {
193ab8d36c9SPeter Brune     *lsnes = fas->next;
194ab8d36c9SPeter Brune     fas = (SNES_FAS *)(*lsnes)->data;
195ab8d36c9SPeter Brune   }
196ab8d36c9SPeter Brune   if (fas->level != level) SETERRQ(((PetscObject)snes)->comm,PETSC_ERR_PLIB,"SNESFAS level hierarchy corrupt");
197ab8d36c9SPeter Brune   PetscFunctionReturn(0);
198ab8d36c9SPeter Brune }
199ab8d36c9SPeter Brune 
200ab8d36c9SPeter Brune #undef __FUNCT__
201ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothUp"
202ab8d36c9SPeter Brune /*@
203ab8d36c9SPeter Brune    SNESFASSetNumberSmoothUp - Sets the number of post-smoothing steps to
204ab8d36c9SPeter Brune    use on all levels.
205ab8d36c9SPeter Brune 
206ab8d36c9SPeter Brune    Logically Collective on SNES
207ab8d36c9SPeter Brune 
208ab8d36c9SPeter Brune    Input Parameters:
209ab8d36c9SPeter Brune +  snes - the multigrid context
210ab8d36c9SPeter Brune -  n    - the number of smoothing steps
211ab8d36c9SPeter Brune 
212ab8d36c9SPeter Brune    Options Database Key:
213ab8d36c9SPeter Brune .  -snes_fas_smoothup <n> - Sets number of pre-smoothing steps
214ab8d36c9SPeter Brune 
215ab8d36c9SPeter Brune    Level: advanced
216ab8d36c9SPeter Brune 
217ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
218ab8d36c9SPeter Brune 
219ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothDown()
220ab8d36c9SPeter Brune @*/
22122d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothUp(SNES snes, PetscInt n)
22222d28d08SBarry Smith {
223ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
22422d28d08SBarry Smith   PetscErrorCode ierr;
22522d28d08SBarry Smith 
226ab8d36c9SPeter Brune   PetscFunctionBegin;
227ab8d36c9SPeter Brune   fas->max_up_it = n;
228656ede7eSPeter Brune   if (!fas->smoothu && fas->level != 0) {
22922d28d08SBarry Smith     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
230ab8d36c9SPeter Brune   }
23122d28d08SBarry Smith   if (fas->smoothu) {
23222d28d08SBarry Smith     ierr = SNESSetTolerances(fas->smoothu, fas->smoothu->abstol, fas->smoothu->rtol, fas->smoothu->stol, n, fas->smoothu->max_funcs);CHKERRQ(ierr);
23322d28d08SBarry Smith   }
234ab8d36c9SPeter Brune   if (fas->next) {
235ab8d36c9SPeter Brune     ierr = SNESFASSetNumberSmoothUp(fas->next, n);CHKERRQ(ierr);
236ab8d36c9SPeter Brune   }
237ab8d36c9SPeter Brune   PetscFunctionReturn(0);
238ab8d36c9SPeter Brune }
239ab8d36c9SPeter Brune 
240ab8d36c9SPeter Brune #undef __FUNCT__
241ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetNumberSmoothDown"
242ab8d36c9SPeter Brune /*@
243ab8d36c9SPeter Brune    SNESFASSetNumberSmoothDown - Sets the number of pre-smoothing steps to
244ab8d36c9SPeter Brune    use on all levels.
245ab8d36c9SPeter Brune 
246ab8d36c9SPeter Brune    Logically Collective on SNES
247ab8d36c9SPeter Brune 
248ab8d36c9SPeter Brune    Input Parameters:
249ab8d36c9SPeter Brune +  snes - the multigrid context
250ab8d36c9SPeter Brune -  n    - the number of smoothing steps
251ab8d36c9SPeter Brune 
252ab8d36c9SPeter Brune    Options Database Key:
253ab8d36c9SPeter Brune .  -snes_fas_smoothdown <n> - Sets number of pre-smoothing steps
254ab8d36c9SPeter Brune 
255ab8d36c9SPeter Brune    Level: advanced
256ab8d36c9SPeter Brune 
257ab8d36c9SPeter Brune .keywords: FAS, MG, smooth, down, pre-smoothing, steps, multigrid
258ab8d36c9SPeter Brune 
259ab8d36c9SPeter Brune .seealso: SNESFASSetNumberSmoothUp()
260ab8d36c9SPeter Brune @*/
26122d28d08SBarry Smith PetscErrorCode SNESFASSetNumberSmoothDown(SNES snes, PetscInt n)
26222d28d08SBarry Smith {
263ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
264ab8d36c9SPeter Brune   PetscErrorCode ierr = 0;
26522d28d08SBarry Smith 
266ab8d36c9SPeter Brune   PetscFunctionBegin;
267ab8d36c9SPeter Brune   if (!fas->smoothd) {
26822d28d08SBarry Smith     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
269ab8d36c9SPeter Brune   }
270ab8d36c9SPeter Brune   ierr = SNESSetTolerances(fas->smoothd, fas->smoothd->abstol, fas->smoothd->rtol, fas->smoothd->stol, n, fas->smoothd->max_funcs);CHKERRQ(ierr);
271ab8d36c9SPeter Brune   fas->max_down_it = n;
272ab8d36c9SPeter Brune   if (fas->next) {
273ab8d36c9SPeter Brune     ierr = SNESFASSetNumberSmoothDown(fas->next, n);CHKERRQ(ierr);
274ab8d36c9SPeter Brune   }
275ab8d36c9SPeter Brune   PetscFunctionReturn(0);
276ab8d36c9SPeter Brune }
277ab8d36c9SPeter Brune 
278ab8d36c9SPeter Brune 
279ab8d36c9SPeter Brune #undef __FUNCT__
280ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetCycles"
281ab8d36c9SPeter Brune /*@
282ab8d36c9SPeter Brune    SNESFASSetCycles - Sets the number of FAS multigrid cycles to use each time a grid is visited.  Use SNESFASSetCyclesOnLevel() for more
283ab8d36c9SPeter Brune    complicated cycling.
284ab8d36c9SPeter Brune 
285ab8d36c9SPeter Brune    Logically Collective on SNES
286ab8d36c9SPeter Brune 
287ab8d36c9SPeter Brune    Input Parameters:
288ab8d36c9SPeter Brune +  snes   - the multigrid context
289ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
290ab8d36c9SPeter Brune 
291ab8d36c9SPeter Brune    Options Database Key:
292ab8d36c9SPeter Brune $  -snes_fas_cycles 1 or 2
293ab8d36c9SPeter Brune 
294ab8d36c9SPeter Brune    Level: advanced
295ab8d36c9SPeter Brune 
296ab8d36c9SPeter Brune .keywords: MG, set, cycles, V-cycle, W-cycle, multigrid
297ab8d36c9SPeter Brune 
298ab8d36c9SPeter Brune .seealso: SNESFASSetCyclesOnLevel()
299ab8d36c9SPeter Brune @*/
30022d28d08SBarry Smith PetscErrorCode SNESFASSetCycles(SNES snes, PetscInt cycles)
30122d28d08SBarry Smith {
302ab8d36c9SPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
303ab8d36c9SPeter Brune   PetscErrorCode ierr;
304ab8d36c9SPeter Brune   PetscBool      isFine;
30522d28d08SBarry Smith 
306ab8d36c9SPeter Brune   PetscFunctionBegin;
30722d28d08SBarry Smith   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
308ab8d36c9SPeter Brune   fas->n_cycles = cycles;
309ab8d36c9SPeter Brune   if (!isFine)
310ab8d36c9SPeter Brune     ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
311ab8d36c9SPeter Brune   if (fas->next) {
312ab8d36c9SPeter Brune     ierr = SNESFASSetCycles(fas->next, cycles);CHKERRQ(ierr);
313ab8d36c9SPeter Brune   }
314ab8d36c9SPeter Brune   PetscFunctionReturn(0);
315ab8d36c9SPeter Brune }
316ab8d36c9SPeter Brune 
317c8c899caSPeter Brune 
318c8c899caSPeter Brune #undef __FUNCT__
319c8c899caSPeter Brune #define __FUNCT__ "SNESFASSetMonitor"
320c8c899caSPeter Brune /*@
321c8c899caSPeter Brune    SNESFASSetMonitor - Sets the method-specific cycle monitoring
322c8c899caSPeter Brune 
323c8c899caSPeter Brune    Logically Collective on SNES
324c8c899caSPeter Brune 
325c8c899caSPeter Brune    Input Parameters:
326c8c899caSPeter Brune +  snes   - the FAS context
327c8c899caSPeter Brune -  flg    - monitor or not
328c8c899caSPeter Brune 
329c8c899caSPeter Brune    Level: advanced
330c8c899caSPeter Brune 
331c8c899caSPeter Brune .keywords: FAS, monitor
332c8c899caSPeter Brune 
333c8c899caSPeter Brune .seealso: SNESFASSetCyclesOnLevel()
334c8c899caSPeter Brune @*/
33522d28d08SBarry Smith PetscErrorCode SNESFASSetMonitor(SNES snes, PetscBool flg)
33622d28d08SBarry Smith {
337c8c899caSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
338c8c899caSPeter Brune   PetscErrorCode ierr;
339c8c899caSPeter Brune   PetscBool      isFine;
340c8c899caSPeter Brune   PetscInt       i, levels = fas->levels;
341c8c899caSPeter Brune   SNES           levelsnes;
34222d28d08SBarry Smith 
343c8c899caSPeter Brune   PetscFunctionBegin;
344c8c899caSPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
345c8c899caSPeter Brune   if (isFine) {
346c8c899caSPeter Brune     for (i = 0; i < levels; i++) {
34722d28d08SBarry Smith       ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
348c8c899caSPeter Brune       fas = (SNES_FAS *)levelsnes->data;
349c8c899caSPeter Brune       if (flg) {
350c8c899caSPeter Brune         fas->monitor = PETSC_VIEWER_STDOUT_(((PetscObject)levelsnes)->comm);CHKERRQ(ierr);
351c8c899caSPeter Brune         /* set the monitors for the upsmoother and downsmoother */
352c8c899caSPeter Brune         ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr);
353c8c899caSPeter Brune         ierr = SNESMonitorSet(levelsnes,SNESMonitorDefault,PETSC_NULL,(PetscErrorCode (*)(void**))PetscViewerDestroy);CHKERRQ(ierr);
354c8c899caSPeter Brune       } else {
355c8c899caSPeter Brune         /* unset the monitors on the coarse levels */
356c8c899caSPeter Brune         if (i != fas->levels - 1) {
357c8c899caSPeter Brune           ierr = SNESMonitorCancel(levelsnes);CHKERRQ(ierr);
358c8c899caSPeter Brune         }
359c8c899caSPeter Brune       }
360c8c899caSPeter Brune     }
361c8c899caSPeter Brune   }
362c8c899caSPeter Brune   PetscFunctionReturn(0);
363c8c899caSPeter Brune }
364c8c899caSPeter Brune 
365ab8d36c9SPeter Brune #undef __FUNCT__
3660dd27c6cSPeter Brune #define __FUNCT__ "SNESFASSetLog"
3670dd27c6cSPeter Brune /*@
3680dd27c6cSPeter Brune    SNESFASSetLog - Sets or unsets time logging for various FAS stages on all levels
3690dd27c6cSPeter Brune 
3700dd27c6cSPeter Brune    Logically Collective on SNES
3710dd27c6cSPeter Brune 
3720dd27c6cSPeter Brune    Input Parameters:
3730dd27c6cSPeter Brune +  snes   - the FAS context
3740dd27c6cSPeter Brune -  flg    - monitor or not
3750dd27c6cSPeter Brune 
3760dd27c6cSPeter Brune    Level: advanced
3770dd27c6cSPeter Brune 
3780dd27c6cSPeter Brune .keywords: FAS, logging
3790dd27c6cSPeter Brune 
3800dd27c6cSPeter Brune .seealso: SNESFASSetMonitor()
3810dd27c6cSPeter Brune @*/
3820dd27c6cSPeter Brune PetscErrorCode SNESFASSetLog(SNES snes, PetscBool flg)
3830dd27c6cSPeter Brune {
3840dd27c6cSPeter Brune   SNES_FAS       *fas = (SNES_FAS *)snes->data;
3850dd27c6cSPeter Brune   PetscErrorCode ierr;
3860dd27c6cSPeter Brune   PetscBool      isFine;
3870dd27c6cSPeter Brune   PetscInt       i, levels = fas->levels;
3880dd27c6cSPeter Brune   SNES           levelsnes;
3890dd27c6cSPeter Brune   char           eventname[128];
3900dd27c6cSPeter Brune 
3910dd27c6cSPeter Brune   PetscFunctionBegin;
3920dd27c6cSPeter Brune   ierr = SNESFASCycleIsFine(snes, &isFine);CHKERRQ(ierr);
3930dd27c6cSPeter Brune   if (isFine) {
3940dd27c6cSPeter Brune     for (i = 0; i < levels; i++) {
3950dd27c6cSPeter Brune       ierr = SNESFASGetCycleSNES(snes, i, &levelsnes);CHKERRQ(ierr);
3960dd27c6cSPeter Brune       fas = (SNES_FAS *)levelsnes->data;
3970dd27c6cSPeter Brune       if (flg) {
3980dd27c6cSPeter Brune         sprintf(eventname,"FASSetup %d",(int)i);
3990dd27c6cSPeter Brune         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsetup);CHKERRQ(ierr);
4000dd27c6cSPeter Brune         sprintf(eventname,"FASSmooth %d",(int)i);
4010dd27c6cSPeter Brune         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventsmoothsolve);CHKERRQ(ierr);
4020dd27c6cSPeter Brune         sprintf(eventname,"FASResid %d",(int)i);
4030dd27c6cSPeter Brune         ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventresidual);CHKERRQ(ierr);
4040dd27c6cSPeter Brune         if (i) {
4050dd27c6cSPeter Brune           sprintf(eventname,"FASInterp %d",(int)i);
4060dd27c6cSPeter Brune           ierr = PetscLogEventRegister(eventname,((PetscObject)snes)->classid,&fas->eventinterprestrict);CHKERRQ(ierr);
4070dd27c6cSPeter Brune         }
4080dd27c6cSPeter Brune       } else {
4090dd27c6cSPeter Brune         fas->eventsmoothsetup    = PETSC_FALSE;
4100dd27c6cSPeter Brune         fas->eventsmoothsolve    = PETSC_FALSE;
4110dd27c6cSPeter Brune         fas->eventresidual       = PETSC_FALSE;
4120dd27c6cSPeter Brune         fas->eventinterprestrict = PETSC_FALSE;
4130dd27c6cSPeter Brune       }
4140dd27c6cSPeter Brune     }
4150dd27c6cSPeter Brune   }
4160dd27c6cSPeter Brune   PetscFunctionReturn(0);
4170dd27c6cSPeter Brune }
4180dd27c6cSPeter Brune 
4190dd27c6cSPeter Brune #undef __FUNCT__
420ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleCreateSmoother_Private"
421ab8d36c9SPeter Brune /*
422ab8d36c9SPeter Brune Creates the default smoother type.
423ab8d36c9SPeter Brune 
42404d7464bSBarry Smith This is SNESNRICHARDSON on each fine level and SNESNEWTONLS on the coarse level.
425ab8d36c9SPeter Brune 
426ab8d36c9SPeter Brune  */
42722d28d08SBarry Smith PetscErrorCode SNESFASCycleCreateSmoother_Private(SNES snes, SNES *smooth)
42822d28d08SBarry Smith {
429ab8d36c9SPeter Brune   SNES_FAS       *fas;
430ab8d36c9SPeter Brune   const char     *optionsprefix;
431ab8d36c9SPeter Brune   char           tprefix[128];
432ab8d36c9SPeter Brune   PetscErrorCode ierr;
433ab8d36c9SPeter Brune   SNES           nsmooth;
43422d28d08SBarry Smith 
435ab8d36c9SPeter Brune   PetscFunctionBegin;
436ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
437ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
438ab8d36c9SPeter Brune   ierr = SNESGetOptionsPrefix(fas->fine, &optionsprefix);CHKERRQ(ierr);
439ab8d36c9SPeter Brune   /* create the default smoother */
440ab8d36c9SPeter Brune   ierr = SNESCreate(((PetscObject)snes)->comm, &nsmooth);CHKERRQ(ierr);
441ab8d36c9SPeter Brune   if (fas->level == 0) {
442ab8d36c9SPeter Brune     sprintf(tprefix,"fas_coarse_");
443ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
444ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
44504d7464bSBarry Smith     ierr = SNESSetType(nsmooth, SNESNEWTONLS);CHKERRQ(ierr);
446e70c42e5SPeter Brune     ierr = SNESSetTolerances(nsmooth, nsmooth->abstol, nsmooth->rtol, nsmooth->stol, nsmooth->max_its, nsmooth->max_funcs);CHKERRQ(ierr);
447ab8d36c9SPeter Brune   } else {
448ab8d36c9SPeter Brune     sprintf(tprefix,"fas_levels_%d_",(int)fas->level);
449ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, optionsprefix);CHKERRQ(ierr);
450ab8d36c9SPeter Brune     ierr = SNESAppendOptionsPrefix(nsmooth, tprefix);CHKERRQ(ierr);
451ab8d36c9SPeter Brune     ierr = SNESSetType(nsmooth, SNESNRICHARDSON);CHKERRQ(ierr);
452e70c42e5SPeter Brune     ierr = SNESSetTolerances(nsmooth, 0.0, 0.0, 0.0, fas->max_down_it, nsmooth->max_funcs);CHKERRQ(ierr);
453ab8d36c9SPeter Brune   }
454ab8d36c9SPeter Brune   ierr = PetscObjectIncrementTabLevel((PetscObject)nsmooth, (PetscObject)snes, 1);CHKERRQ(ierr);
455f89ba88eSPeter Brune   ierr = PetscLogObjectParent(snes,nsmooth);CHKERRQ(ierr);
456f89ba88eSPeter Brune   ierr = PetscObjectCopyFortranFunctionPointers((PetscObject)snes, (PetscObject)nsmooth);CHKERRQ(ierr);
457ab8d36c9SPeter Brune   *smooth = nsmooth;
458ab8d36c9SPeter Brune   PetscFunctionReturn(0);
459ab8d36c9SPeter Brune }
460ab8d36c9SPeter Brune 
461ab8d36c9SPeter Brune /* ------------- Functions called on a particular level ----------------- */
462ab8d36c9SPeter Brune 
463ab8d36c9SPeter Brune #undef __FUNCT__
464ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleSetCycles"
465ab8d36c9SPeter Brune /*@
466ab8d36c9SPeter Brune    SNESFASCycleSetCycles - Sets the number of cycles on a particular level.
467ab8d36c9SPeter Brune 
468ab8d36c9SPeter Brune    Logically Collective on SNES
469ab8d36c9SPeter Brune 
470ab8d36c9SPeter Brune    Input Parameters:
471ab8d36c9SPeter Brune +  snes   - the multigrid context
472ab8d36c9SPeter Brune .  level  - the level to set the number of cycles on
473ab8d36c9SPeter Brune -  cycles - the number of cycles -- 1 for V-cycle, 2 for W-cycle
474ab8d36c9SPeter Brune 
475ab8d36c9SPeter Brune    Level: advanced
476ab8d36c9SPeter Brune 
477ab8d36c9SPeter Brune .keywords: SNES, FAS, set, cycles, V-cycle, W-cycle, multigrid
478ab8d36c9SPeter Brune 
479ab8d36c9SPeter Brune .seealso: SNESFASSetCycles()
480ab8d36c9SPeter Brune @*/
48122d28d08SBarry Smith PetscErrorCode SNESFASCycleSetCycles(SNES snes, PetscInt cycles)
48222d28d08SBarry Smith {
483ab8d36c9SPeter Brune   SNES_FAS       *fas =  (SNES_FAS *)snes->data;
484ab8d36c9SPeter Brune   PetscErrorCode ierr;
48522d28d08SBarry Smith 
486ab8d36c9SPeter Brune   PetscFunctionBegin;
487ab8d36c9SPeter Brune   fas->n_cycles = cycles;
488ab8d36c9SPeter Brune   ierr = SNESSetTolerances(snes, snes->abstol, snes->rtol, snes->stol, cycles, snes->max_funcs);CHKERRQ(ierr);
489ab8d36c9SPeter Brune   PetscFunctionReturn(0);
490ab8d36c9SPeter Brune }
491ab8d36c9SPeter Brune 
492ab8d36c9SPeter Brune 
493ab8d36c9SPeter Brune #undef __FUNCT__
494ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmoother"
495ab8d36c9SPeter Brune /*@
496ab8d36c9SPeter Brune    SNESFASCycleGetSmoother - Gets the smoother on a particular cycle level.
497ab8d36c9SPeter Brune 
498ab8d36c9SPeter Brune    Logically Collective on SNES
499ab8d36c9SPeter Brune 
500ab8d36c9SPeter Brune    Input Parameters:
501ab8d36c9SPeter Brune .  snes   - the multigrid context
502ab8d36c9SPeter Brune 
503ab8d36c9SPeter Brune    Output Parameters:
504ab8d36c9SPeter Brune .  smooth - the smoother
505ab8d36c9SPeter Brune 
506ab8d36c9SPeter Brune    Level: advanced
507ab8d36c9SPeter Brune 
508ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
509ab8d36c9SPeter Brune 
510ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmootherDown()
511ab8d36c9SPeter Brune @*/
512ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmoother(SNES snes, SNES *smooth)
513ab8d36c9SPeter Brune {
514ab8d36c9SPeter Brune   SNES_FAS       *fas;
515*5fd66863SKarl Rupp 
516ab8d36c9SPeter Brune   PetscFunctionBegin;
517ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
518ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
519ab8d36c9SPeter Brune   *smooth = fas->smoothd;
520ab8d36c9SPeter Brune   PetscFunctionReturn(0);
521ab8d36c9SPeter Brune }
522ab8d36c9SPeter Brune #undef __FUNCT__
523ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherUp"
524ab8d36c9SPeter Brune /*@
525ab8d36c9SPeter Brune    SNESFASCycleGetSmootherUp - Gets the up smoother on a particular cycle level.
526ab8d36c9SPeter Brune 
527ab8d36c9SPeter Brune    Logically Collective on SNES
528ab8d36c9SPeter Brune 
529ab8d36c9SPeter Brune    Input Parameters:
530ab8d36c9SPeter Brune .  snes   - the multigrid context
531ab8d36c9SPeter Brune 
532ab8d36c9SPeter Brune    Output Parameters:
533ab8d36c9SPeter Brune .  smoothu - the smoother
534ab8d36c9SPeter Brune 
535ab8d36c9SPeter Brune    Notes:
536ab8d36c9SPeter Brune    Returns the downsmoother if no up smoother is available.  This enables transparent
537ab8d36c9SPeter Brune    default behavior in the process of the solve.
538ab8d36c9SPeter Brune 
539ab8d36c9SPeter Brune    Level: advanced
540ab8d36c9SPeter Brune 
541ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
542ab8d36c9SPeter Brune 
543ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmoother(), SNESFASCycleGetSmootherDown()
544ab8d36c9SPeter Brune @*/
545ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherUp(SNES snes, SNES *smoothu)
546ab8d36c9SPeter Brune {
547ab8d36c9SPeter Brune   SNES_FAS       *fas;
548*5fd66863SKarl Rupp 
549ab8d36c9SPeter Brune   PetscFunctionBegin;
550ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
551ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
552ab8d36c9SPeter Brune   if (!fas->smoothu) {
553ab8d36c9SPeter Brune     *smoothu = fas->smoothd;
554ab8d36c9SPeter Brune   } else {
555ab8d36c9SPeter Brune     *smoothu = fas->smoothu;
556ab8d36c9SPeter Brune   }
557ab8d36c9SPeter Brune   PetscFunctionReturn(0);
558ab8d36c9SPeter Brune }
559ab8d36c9SPeter Brune 
560ab8d36c9SPeter Brune #undef __FUNCT__
561ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetSmootherDown"
562ab8d36c9SPeter Brune /*@
563ab8d36c9SPeter Brune    SNESFASCycleGetSmootherDown - Gets the down smoother on a particular cycle level.
564ab8d36c9SPeter Brune 
565ab8d36c9SPeter Brune    Logically Collective on SNES
566ab8d36c9SPeter Brune 
567ab8d36c9SPeter Brune    Input Parameters:
568ab8d36c9SPeter Brune .  snes   - the multigrid context
569ab8d36c9SPeter Brune 
570ab8d36c9SPeter Brune    Output Parameters:
571ab8d36c9SPeter Brune .  smoothd - the smoother
572ab8d36c9SPeter Brune 
573ab8d36c9SPeter Brune    Level: advanced
574ab8d36c9SPeter Brune 
575ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
576ab8d36c9SPeter Brune 
577ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
578ab8d36c9SPeter Brune @*/
579ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetSmootherDown(SNES snes, SNES *smoothd)
580ab8d36c9SPeter Brune {
581ab8d36c9SPeter Brune   SNES_FAS       *fas;
582*5fd66863SKarl Rupp 
583ab8d36c9SPeter Brune   PetscFunctionBegin;
584ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
585ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
586ab8d36c9SPeter Brune   *smoothd = fas->smoothd;
587ab8d36c9SPeter Brune   PetscFunctionReturn(0);
588ab8d36c9SPeter Brune }
589ab8d36c9SPeter Brune 
590ab8d36c9SPeter Brune 
591ab8d36c9SPeter Brune #undef __FUNCT__
592ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetCorrection"
593ab8d36c9SPeter Brune /*@
594ab8d36c9SPeter Brune    SNESFASCycleGetCorrection - Gets the coarse correction FAS context for this level
595ab8d36c9SPeter Brune 
596ab8d36c9SPeter Brune    Logically Collective on SNES
597ab8d36c9SPeter Brune 
598ab8d36c9SPeter Brune    Input Parameters:
599ab8d36c9SPeter Brune .  snes   - the multigrid context
600ab8d36c9SPeter Brune 
601ab8d36c9SPeter Brune    Output Parameters:
602ab8d36c9SPeter Brune .  correction - the coarse correction on this level
603ab8d36c9SPeter Brune 
604ab8d36c9SPeter Brune    Notes:
605ab8d36c9SPeter Brune    Returns PETSC_NULL on the coarsest level.
606ab8d36c9SPeter Brune 
607ab8d36c9SPeter Brune    Level: advanced
608ab8d36c9SPeter Brune 
609ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
610ab8d36c9SPeter Brune 
611ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
612ab8d36c9SPeter Brune @*/
613ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetCorrection(SNES snes, SNES *correction)
614ab8d36c9SPeter Brune {
615ab8d36c9SPeter Brune   SNES_FAS       *fas;
616*5fd66863SKarl Rupp 
617ab8d36c9SPeter Brune   PetscFunctionBegin;
618ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
619ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
620ab8d36c9SPeter Brune   *correction = fas->next;
621ab8d36c9SPeter Brune   PetscFunctionReturn(0);
622ab8d36c9SPeter Brune }
623ab8d36c9SPeter Brune 
624ab8d36c9SPeter Brune #undef __FUNCT__
625ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInterpolation"
626ab8d36c9SPeter Brune /*@
627ab8d36c9SPeter Brune    SNESFASCycleGetInterpolation - Gets the interpolation on this level
628ab8d36c9SPeter Brune 
629ab8d36c9SPeter Brune    Logically Collective on SNES
630ab8d36c9SPeter Brune 
631ab8d36c9SPeter Brune    Input Parameters:
632ab8d36c9SPeter Brune .  snes   - the multigrid context
633ab8d36c9SPeter Brune 
634ab8d36c9SPeter Brune    Output Parameters:
635ab8d36c9SPeter Brune .  mat    - the interpolation operator on this level
636ab8d36c9SPeter Brune 
637ab8d36c9SPeter Brune    Level: developer
638ab8d36c9SPeter Brune 
639ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
640ab8d36c9SPeter Brune 
641ab8d36c9SPeter Brune .seealso: SNESFASCycleGetSmootherUp(), SNESFASCycleGetSmoother()
642ab8d36c9SPeter Brune @*/
643ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInterpolation(SNES snes, Mat *mat)
644ab8d36c9SPeter Brune {
645ab8d36c9SPeter Brune   SNES_FAS       *fas;
646*5fd66863SKarl Rupp 
647ab8d36c9SPeter Brune   PetscFunctionBegin;
648ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
649ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
650ab8d36c9SPeter Brune   *mat = fas->interpolate;
651ab8d36c9SPeter Brune   PetscFunctionReturn(0);
652ab8d36c9SPeter Brune }
653ab8d36c9SPeter Brune 
654ab8d36c9SPeter Brune 
655ab8d36c9SPeter Brune #undef __FUNCT__
656ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRestriction"
657ab8d36c9SPeter Brune /*@
658ab8d36c9SPeter Brune    SNESFASCycleGetRestriction - Gets the restriction on this level
659ab8d36c9SPeter Brune 
660ab8d36c9SPeter Brune    Logically Collective on SNES
661ab8d36c9SPeter Brune 
662ab8d36c9SPeter Brune    Input Parameters:
663ab8d36c9SPeter Brune .  snes   - the multigrid context
664ab8d36c9SPeter Brune 
665ab8d36c9SPeter Brune    Output Parameters:
666ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
667ab8d36c9SPeter Brune 
668ab8d36c9SPeter Brune    Level: developer
669ab8d36c9SPeter Brune 
670ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
671ab8d36c9SPeter Brune 
672ab8d36c9SPeter Brune .seealso: SNESFASGetRestriction(), SNESFASCycleGetInterpolation()
673ab8d36c9SPeter Brune @*/
674ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRestriction(SNES snes, Mat *mat)
675ab8d36c9SPeter Brune {
676ab8d36c9SPeter Brune   SNES_FAS       *fas;
677*5fd66863SKarl Rupp 
678ab8d36c9SPeter Brune   PetscFunctionBegin;
679ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
680ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
681ab8d36c9SPeter Brune   *mat = fas->restrct;
682ab8d36c9SPeter Brune   PetscFunctionReturn(0);
683ab8d36c9SPeter Brune }
684ab8d36c9SPeter Brune 
685ab8d36c9SPeter Brune 
686ab8d36c9SPeter Brune #undef __FUNCT__
687ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetInjection"
688ab8d36c9SPeter Brune /*@
689ab8d36c9SPeter Brune    SNESFASCycleGetInjection - Gets the injection on this level
690ab8d36c9SPeter Brune 
691ab8d36c9SPeter Brune    Logically Collective on SNES
692ab8d36c9SPeter Brune 
693ab8d36c9SPeter Brune    Input Parameters:
694ab8d36c9SPeter Brune .  snes   - the multigrid context
695ab8d36c9SPeter Brune 
696ab8d36c9SPeter Brune    Output Parameters:
697ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
698ab8d36c9SPeter Brune 
699ab8d36c9SPeter Brune    Level: developer
700ab8d36c9SPeter Brune 
701ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
702ab8d36c9SPeter Brune 
703ab8d36c9SPeter Brune .seealso: SNESFASGetInjection(), SNESFASCycleGetRestriction()
704ab8d36c9SPeter Brune @*/
705ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetInjection(SNES snes, Mat *mat)
706ab8d36c9SPeter Brune {
707ab8d36c9SPeter Brune   SNES_FAS       *fas;
708*5fd66863SKarl Rupp 
709ab8d36c9SPeter Brune   PetscFunctionBegin;
710ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
711ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
712ab8d36c9SPeter Brune   *mat = fas->inject;
713ab8d36c9SPeter Brune   PetscFunctionReturn(0);
714ab8d36c9SPeter Brune }
715ab8d36c9SPeter Brune 
716ab8d36c9SPeter Brune #undef __FUNCT__
717ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleGetRScale"
718ab8d36c9SPeter Brune /*@
719ab8d36c9SPeter Brune    SNESFASCycleGetRScale - Gets the injection on this level
720ab8d36c9SPeter Brune 
721ab8d36c9SPeter Brune    Logically Collective on SNES
722ab8d36c9SPeter Brune 
723ab8d36c9SPeter Brune    Input Parameters:
724ab8d36c9SPeter Brune .  snes   - the multigrid context
725ab8d36c9SPeter Brune 
726ab8d36c9SPeter Brune    Output Parameters:
727ab8d36c9SPeter Brune .  mat    - the restriction operator on this level
728ab8d36c9SPeter Brune 
729ab8d36c9SPeter Brune    Level: developer
730ab8d36c9SPeter Brune 
731ab8d36c9SPeter Brune .keywords: SNES, FAS, get, smoother, multigrid
732ab8d36c9SPeter Brune 
733ab8d36c9SPeter Brune .seealso: SNESFASCycleGetRestriction(), SNESFASGetRScale()
734ab8d36c9SPeter Brune @*/
735ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleGetRScale(SNES snes, Vec *vec)
736ab8d36c9SPeter Brune {
737ab8d36c9SPeter Brune   SNES_FAS       *fas;
738*5fd66863SKarl Rupp 
739ab8d36c9SPeter Brune   PetscFunctionBegin;
740ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
741ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
742ab8d36c9SPeter Brune   *vec = fas->rscale;
743ab8d36c9SPeter Brune   PetscFunctionReturn(0);
744ab8d36c9SPeter Brune }
745ab8d36c9SPeter Brune 
746ab8d36c9SPeter Brune #undef __FUNCT__
747ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASCycleIsFine"
748ab8d36c9SPeter Brune /*@
749ab8d36c9SPeter Brune    SNESFASCycleIsFine - Determines if a given cycle is the fine level.
750ab8d36c9SPeter Brune 
751ab8d36c9SPeter Brune    Logically Collective on SNES
752ab8d36c9SPeter Brune 
753ab8d36c9SPeter Brune    Input Parameters:
754ab8d36c9SPeter Brune .  snes   - the FAS context
755ab8d36c9SPeter Brune 
756ab8d36c9SPeter Brune    Output Parameters:
757ab8d36c9SPeter Brune .  flg - indicates if this is the fine level or not
758ab8d36c9SPeter Brune 
759ab8d36c9SPeter Brune    Level: advanced
760ab8d36c9SPeter Brune 
761ab8d36c9SPeter Brune .keywords: SNES, FAS
762ab8d36c9SPeter Brune 
763ab8d36c9SPeter Brune .seealso: SNESFASSetLevels()
764ab8d36c9SPeter Brune @*/
765ab8d36c9SPeter Brune PetscErrorCode SNESFASCycleIsFine(SNES snes, PetscBool *flg)
766ab8d36c9SPeter Brune {
767ab8d36c9SPeter Brune   SNES_FAS       *fas;
768*5fd66863SKarl Rupp 
769ab8d36c9SPeter Brune   PetscFunctionBegin;
770ab8d36c9SPeter Brune   PetscValidHeaderSpecific(snes,SNES_CLASSID,1);
771ab8d36c9SPeter Brune   fas = (SNES_FAS*)snes->data;
772ab8d36c9SPeter Brune   if (fas->level == fas->levels - 1) {
773ab8d36c9SPeter Brune     *flg = PETSC_TRUE;
774ab8d36c9SPeter Brune   } else {
775ab8d36c9SPeter Brune     *flg = PETSC_FALSE;
776ab8d36c9SPeter Brune   }
777ab8d36c9SPeter Brune   PetscFunctionReturn(0);
778ab8d36c9SPeter Brune }
779ab8d36c9SPeter Brune 
780ab8d36c9SPeter Brune /* ---------- functions called on the finest level that return level-specific information ---------- */
781ab8d36c9SPeter Brune 
782ab8d36c9SPeter Brune #undef __FUNCT__
783ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInterpolation"
784ab8d36c9SPeter Brune /*@
785ab8d36c9SPeter Brune    SNESFASSetInterpolation - Sets the function to be used to calculate the
786ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
787ab8d36c9SPeter Brune 
788ab8d36c9SPeter Brune    Input Parameters:
789ab8d36c9SPeter Brune +  snes      - the multigrid context
790ab8d36c9SPeter Brune .  mat       - the interpolation operator
791ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
792ab8d36c9SPeter Brune 
793ab8d36c9SPeter Brune    Level: advanced
794ab8d36c9SPeter Brune 
795ab8d36c9SPeter Brune    Notes:
796ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the restriction
797ab8d36c9SPeter Brune     for the same level.
798ab8d36c9SPeter Brune 
799ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
800ab8d36c9SPeter Brune     out from the matrix size which one it is.
801ab8d36c9SPeter Brune 
802ab8d36c9SPeter Brune .keywords:  FAS, multigrid, set, interpolate, level
803ab8d36c9SPeter Brune 
804ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction(), SNESFASSetRScale()
805ab8d36c9SPeter Brune @*/
8060adebc6cSBarry Smith PetscErrorCode SNESFASSetInterpolation(SNES snes, PetscInt level, Mat mat)
8070adebc6cSBarry Smith {
80822d28d08SBarry Smith   SNES_FAS       *fas;
809ab8d36c9SPeter Brune   PetscErrorCode ierr;
810ab8d36c9SPeter Brune   SNES           levelsnes;
81122d28d08SBarry Smith 
812ab8d36c9SPeter Brune   PetscFunctionBegin;
813ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
814ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
815ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
816ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->interpolate);CHKERRQ(ierr);
817ab8d36c9SPeter Brune   fas->interpolate = mat;
818ab8d36c9SPeter Brune   PetscFunctionReturn(0);
819ab8d36c9SPeter Brune }
820ab8d36c9SPeter Brune 
821ab8d36c9SPeter Brune #undef __FUNCT__
822ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInterpolation"
823ab8d36c9SPeter Brune /*@
824ab8d36c9SPeter Brune    SNESFASGetInterpolation - Gets the matrix used to calculate the
825ab8d36c9SPeter Brune    interpolation from l-1 to the lth level
826ab8d36c9SPeter Brune 
827ab8d36c9SPeter Brune    Input Parameters:
828ab8d36c9SPeter Brune +  snes      - the multigrid context
829ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
830ab8d36c9SPeter Brune 
831ab8d36c9SPeter Brune    Output Parameters:
832ab8d36c9SPeter Brune .  mat       - the interpolation operator
833ab8d36c9SPeter Brune 
834ab8d36c9SPeter Brune    Level: advanced
835ab8d36c9SPeter Brune 
836ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, interpolate, level
837ab8d36c9SPeter Brune 
838ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASGetInjection(), SNESFASGetRestriction(), SNESFASGetRScale()
839ab8d36c9SPeter Brune @*/
84022d28d08SBarry Smith PetscErrorCode SNESFASGetInterpolation(SNES snes, PetscInt level, Mat *mat)
84122d28d08SBarry Smith {
84222d28d08SBarry Smith   SNES_FAS       *fas;
843ab8d36c9SPeter Brune   PetscErrorCode ierr;
844ab8d36c9SPeter Brune   SNES           levelsnes;
84522d28d08SBarry Smith 
846ab8d36c9SPeter Brune   PetscFunctionBegin;
847ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
848ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
849ab8d36c9SPeter Brune   *mat = fas->interpolate;
850ab8d36c9SPeter Brune   PetscFunctionReturn(0);
851ab8d36c9SPeter Brune }
852ab8d36c9SPeter Brune 
853ab8d36c9SPeter Brune #undef __FUNCT__
854ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRestriction"
855ab8d36c9SPeter Brune /*@
856ab8d36c9SPeter Brune    SNESFASSetRestriction - Sets the function to be used to restrict the defect
857ab8d36c9SPeter Brune    from level l to l-1.
858ab8d36c9SPeter Brune 
859ab8d36c9SPeter Brune    Input Parameters:
860ab8d36c9SPeter Brune +  snes  - the multigrid context
861ab8d36c9SPeter Brune .  mat   - the restriction matrix
862ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
863ab8d36c9SPeter Brune 
864ab8d36c9SPeter Brune    Level: advanced
865ab8d36c9SPeter Brune 
866ab8d36c9SPeter Brune    Notes:
867ab8d36c9SPeter Brune           Usually this is the same matrix used also to set the interpolation
868ab8d36c9SPeter Brune     for the same level.
869ab8d36c9SPeter Brune 
870ab8d36c9SPeter Brune           One can pass in the interpolation matrix or its transpose; PETSc figures
871ab8d36c9SPeter Brune     out from the matrix size which one it is.
872ab8d36c9SPeter Brune 
873ab8d36c9SPeter Brune          If you do not set this, the transpose of the Mat set with SNESFASSetInterpolation()
874ab8d36c9SPeter Brune     is used.
875ab8d36c9SPeter Brune 
876ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
877ab8d36c9SPeter Brune 
878ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetInjection()
879ab8d36c9SPeter Brune @*/
88022d28d08SBarry Smith PetscErrorCode SNESFASSetRestriction(SNES snes, PetscInt level, Mat mat)
88122d28d08SBarry Smith {
88222d28d08SBarry Smith   SNES_FAS       *fas;
883ab8d36c9SPeter Brune   PetscErrorCode ierr;
884ab8d36c9SPeter Brune   SNES           levelsnes;
88522d28d08SBarry Smith 
886ab8d36c9SPeter Brune   PetscFunctionBegin;
887ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
888ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
889ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
890ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->restrct);CHKERRQ(ierr);
891ab8d36c9SPeter Brune   fas->restrct = mat;
892ab8d36c9SPeter Brune   PetscFunctionReturn(0);
893ab8d36c9SPeter Brune }
894ab8d36c9SPeter Brune 
895ab8d36c9SPeter Brune #undef __FUNCT__
896ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetRestriction"
897ab8d36c9SPeter Brune /*@
898ab8d36c9SPeter Brune    SNESFASGetRestriction - Gets the matrix used to calculate the
899ab8d36c9SPeter Brune    restriction from l to the l-1th level
900ab8d36c9SPeter Brune 
901ab8d36c9SPeter Brune    Input Parameters:
902ab8d36c9SPeter Brune +  snes      - the multigrid context
903ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
904ab8d36c9SPeter Brune 
905ab8d36c9SPeter Brune    Output Parameters:
906ab8d36c9SPeter Brune .  mat       - the interpolation operator
907ab8d36c9SPeter Brune 
908ab8d36c9SPeter Brune    Level: advanced
909ab8d36c9SPeter Brune 
910ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, restrict, level
911ab8d36c9SPeter Brune 
912ab8d36c9SPeter Brune .seealso: SNESFASSetRestriction(), SNESFASGetInjection(), SNESFASGetInterpolation(), SNESFASGetRScale()
913ab8d36c9SPeter Brune @*/
91422d28d08SBarry Smith PetscErrorCode SNESFASGetRestriction(SNES snes, PetscInt level, Mat *mat)
91522d28d08SBarry Smith {
91622d28d08SBarry Smith   SNES_FAS       *fas;
917ab8d36c9SPeter Brune   PetscErrorCode ierr;
918ab8d36c9SPeter Brune   SNES           levelsnes;
91922d28d08SBarry Smith 
920ab8d36c9SPeter Brune   PetscFunctionBegin;
921ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
922ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
923ab8d36c9SPeter Brune   *mat = fas->restrct;
924ab8d36c9SPeter Brune   PetscFunctionReturn(0);
925ab8d36c9SPeter Brune }
926ab8d36c9SPeter Brune 
927ab8d36c9SPeter Brune 
928ab8d36c9SPeter Brune #undef __FUNCT__
929ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetInjection"
930ab8d36c9SPeter Brune /*@
931ab8d36c9SPeter Brune    SNESFASSetInjection - Sets the function to be used to inject the solution
932ab8d36c9SPeter Brune    from level l to l-1.
933ab8d36c9SPeter Brune 
934ab8d36c9SPeter Brune    Input Parameters:
935ab8d36c9SPeter Brune  +  snes  - the multigrid context
936ab8d36c9SPeter Brune .  mat   - the restriction matrix
937ab8d36c9SPeter Brune -  level - the level (0 is coarsest) to supply [Do not supply 0]
938ab8d36c9SPeter Brune 
939ab8d36c9SPeter Brune    Level: advanced
940ab8d36c9SPeter Brune 
941ab8d36c9SPeter Brune    Notes:
942ab8d36c9SPeter Brune          If you do not set this, the restriction and rscale is used to
943ab8d36c9SPeter Brune    project the solution instead.
944ab8d36c9SPeter Brune 
945ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
946ab8d36c9SPeter Brune 
947ab8d36c9SPeter Brune .seealso: SNESFASSetInterpolation(), SNESFASSetRestriction()
948ab8d36c9SPeter Brune @*/
94922d28d08SBarry Smith PetscErrorCode SNESFASSetInjection(SNES snes, PetscInt level, Mat mat)
95022d28d08SBarry Smith {
95122d28d08SBarry Smith   SNES_FAS       *fas;
952ab8d36c9SPeter Brune   PetscErrorCode ierr;
953ab8d36c9SPeter Brune   SNES           levelsnes;
95422d28d08SBarry Smith 
955ab8d36c9SPeter Brune   PetscFunctionBegin;
956ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
957ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
958ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)mat);CHKERRQ(ierr);
959ab8d36c9SPeter Brune   ierr = MatDestroy(&fas->inject);CHKERRQ(ierr);
960ab8d36c9SPeter Brune   fas->inject = mat;
961ab8d36c9SPeter Brune   PetscFunctionReturn(0);
962ab8d36c9SPeter Brune }
963ab8d36c9SPeter Brune 
964ab8d36c9SPeter Brune 
965ab8d36c9SPeter Brune #undef __FUNCT__
966ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetInjection"
967ab8d36c9SPeter Brune /*@
968ab8d36c9SPeter Brune    SNESFASGetInjection - Gets the matrix used to calculate the
969ab8d36c9SPeter Brune    injection from l-1 to the lth level
970ab8d36c9SPeter Brune 
971ab8d36c9SPeter Brune    Input Parameters:
972ab8d36c9SPeter Brune +  snes      - the multigrid context
973ab8d36c9SPeter Brune -  level     - the level (0 is coarsest) to supply [do not supply 0]
974ab8d36c9SPeter Brune 
975ab8d36c9SPeter Brune    Output Parameters:
976ab8d36c9SPeter Brune .  mat       - the injection operator
977ab8d36c9SPeter Brune 
978ab8d36c9SPeter Brune    Level: advanced
979ab8d36c9SPeter Brune 
980ab8d36c9SPeter Brune .keywords:  FAS, multigrid, get, restrict, level
981ab8d36c9SPeter Brune 
982ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASGetRestriction(), SNESFASGetInterpolation(), SNESFASGetRScale()
983ab8d36c9SPeter Brune @*/
98422d28d08SBarry Smith PetscErrorCode SNESFASGetInjection(SNES snes, PetscInt level, Mat *mat)
98522d28d08SBarry Smith {
98622d28d08SBarry Smith   SNES_FAS       *fas;
987ab8d36c9SPeter Brune   PetscErrorCode ierr;
988ab8d36c9SPeter Brune   SNES           levelsnes;
98922d28d08SBarry Smith 
990ab8d36c9SPeter Brune   PetscFunctionBegin;
991ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
992ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
993ab8d36c9SPeter Brune   *mat = fas->inject;
994ab8d36c9SPeter Brune   PetscFunctionReturn(0);
995ab8d36c9SPeter Brune }
996ab8d36c9SPeter Brune 
997ab8d36c9SPeter Brune #undef __FUNCT__
998ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASSetRScale"
999ab8d36c9SPeter Brune /*@
1000ab8d36c9SPeter Brune    SNESFASSetRScale - Sets the scaling factor of the restriction
1001ab8d36c9SPeter Brune    operator from level l to l-1.
1002ab8d36c9SPeter Brune 
1003ab8d36c9SPeter Brune    Input Parameters:
1004ab8d36c9SPeter Brune +  snes   - the multigrid context
1005ab8d36c9SPeter Brune .  rscale - the restriction scaling
1006ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply [Do not supply 0]
1007ab8d36c9SPeter Brune 
1008ab8d36c9SPeter Brune    Level: advanced
1009ab8d36c9SPeter Brune 
1010ab8d36c9SPeter Brune    Notes:
1011ab8d36c9SPeter Brune          This is only used in the case that the injection is not set.
1012ab8d36c9SPeter Brune 
1013ab8d36c9SPeter Brune .keywords: FAS, MG, set, multigrid, restriction, level
1014ab8d36c9SPeter Brune 
1015ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1016ab8d36c9SPeter Brune @*/
101722d28d08SBarry Smith PetscErrorCode SNESFASSetRScale(SNES snes, PetscInt level, Vec rscale)
101822d28d08SBarry Smith {
1019ab8d36c9SPeter Brune   SNES_FAS       *fas;
1020ab8d36c9SPeter Brune   PetscErrorCode ierr;
1021ab8d36c9SPeter Brune   SNES           levelsnes;
102222d28d08SBarry Smith 
1023ab8d36c9SPeter Brune   PetscFunctionBegin;
1024ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1025ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1026ab8d36c9SPeter Brune   ierr = PetscObjectReference((PetscObject)rscale);CHKERRQ(ierr);
1027ab8d36c9SPeter Brune   ierr = VecDestroy(&fas->rscale);CHKERRQ(ierr);
1028ab8d36c9SPeter Brune   fas->rscale = rscale;
1029ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1030ab8d36c9SPeter Brune }
1031ab8d36c9SPeter Brune 
1032ab8d36c9SPeter Brune #undef __FUNCT__
1033ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmoother"
1034ab8d36c9SPeter Brune /*@
1035ab8d36c9SPeter Brune    SNESFASGetSmoother - Gets the default smoother on a level.
1036ab8d36c9SPeter Brune 
1037ab8d36c9SPeter Brune    Input Parameters:
1038ab8d36c9SPeter Brune +  snes   - the multigrid context
1039ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
1040ab8d36c9SPeter Brune 
1041ab8d36c9SPeter Brune    Output Parameters:
1042ab8d36c9SPeter Brune    smooth  - the smoother
1043ab8d36c9SPeter Brune 
1044ab8d36c9SPeter Brune    Level: advanced
1045ab8d36c9SPeter Brune 
1046ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
1047ab8d36c9SPeter Brune 
1048ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1049ab8d36c9SPeter Brune @*/
105022d28d08SBarry Smith PetscErrorCode SNESFASGetSmoother(SNES snes, PetscInt level, SNES *smooth)
105122d28d08SBarry Smith {
1052ab8d36c9SPeter Brune   SNES_FAS       *fas;
1053ab8d36c9SPeter Brune   PetscErrorCode ierr;
1054ab8d36c9SPeter Brune   SNES           levelsnes;
105522d28d08SBarry Smith 
1056ab8d36c9SPeter Brune   PetscFunctionBegin;
1057ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1058ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1059ab8d36c9SPeter Brune   if (!fas->smoothd) {
1060ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1061ab8d36c9SPeter Brune   }
1062ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1063ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1064ab8d36c9SPeter Brune }
1065ab8d36c9SPeter Brune 
1066ab8d36c9SPeter Brune #undef __FUNCT__
1067ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherDown"
1068ab8d36c9SPeter Brune /*@
1069ab8d36c9SPeter Brune    SNESFASGetSmootherDown - Gets the downsmoother on a level.
1070ab8d36c9SPeter Brune 
1071ab8d36c9SPeter Brune    Input Parameters:
1072ab8d36c9SPeter Brune +  snes   - the multigrid context
1073ab8d36c9SPeter Brune -  level  - the level (0 is coarsest) to supply
1074ab8d36c9SPeter Brune 
1075ab8d36c9SPeter Brune    Output Parameters:
1076ab8d36c9SPeter Brune    smooth  - the smoother
1077ab8d36c9SPeter Brune 
1078ab8d36c9SPeter Brune    Level: advanced
1079ab8d36c9SPeter Brune 
1080ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
1081ab8d36c9SPeter Brune 
1082ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1083ab8d36c9SPeter Brune @*/
108422d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherDown(SNES snes, PetscInt level, SNES *smooth)
108522d28d08SBarry Smith {
1086ab8d36c9SPeter Brune   SNES_FAS       *fas;
1087ab8d36c9SPeter Brune   PetscErrorCode ierr;
1088ab8d36c9SPeter Brune   SNES           levelsnes;
108922d28d08SBarry Smith 
1090ab8d36c9SPeter Brune   PetscFunctionBegin;
1091ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1092ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1093ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1094ab8d36c9SPeter Brune   if (!fas->smoothd) {
1095ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1096ab8d36c9SPeter Brune   }
1097ab8d36c9SPeter Brune   if (!fas->smoothu) {
1098ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
1099ab8d36c9SPeter Brune   }
1100ab8d36c9SPeter Brune   *smooth = fas->smoothd;
1101ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1102ab8d36c9SPeter Brune }
1103ab8d36c9SPeter Brune 
1104ab8d36c9SPeter Brune #undef __FUNCT__
1105ab8d36c9SPeter Brune #define __FUNCT__ "SNESFASGetSmootherUp"
1106ab8d36c9SPeter Brune /*@
1107ab8d36c9SPeter Brune    SNESFASGetSmootherUp - Gets the upsmoother on a level.
1108ab8d36c9SPeter Brune 
1109ab8d36c9SPeter Brune    Input Parameters:
1110ab8d36c9SPeter Brune +  snes   - the multigrid context
1111ab8d36c9SPeter Brune -  level  - the level (0 is coarsest)
1112ab8d36c9SPeter Brune 
1113ab8d36c9SPeter Brune    Output Parameters:
1114ab8d36c9SPeter Brune    smooth  - the smoother
1115ab8d36c9SPeter Brune 
1116ab8d36c9SPeter Brune    Level: advanced
1117ab8d36c9SPeter Brune 
1118ab8d36c9SPeter Brune .keywords: FAS, MG, get, multigrid, smoother, level
1119ab8d36c9SPeter Brune 
1120ab8d36c9SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1121ab8d36c9SPeter Brune @*/
112222d28d08SBarry Smith PetscErrorCode SNESFASGetSmootherUp(SNES snes, PetscInt level, SNES *smooth)
112322d28d08SBarry Smith {
1124ab8d36c9SPeter Brune   SNES_FAS       *fas;
1125ab8d36c9SPeter Brune   PetscErrorCode ierr;
1126ab8d36c9SPeter Brune   SNES           levelsnes;
112722d28d08SBarry Smith 
1128ab8d36c9SPeter Brune   PetscFunctionBegin;
1129ab8d36c9SPeter Brune   ierr = SNESFASGetCycleSNES(snes, level, &levelsnes);CHKERRQ(ierr);
1130ab8d36c9SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1131ab8d36c9SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1132ab8d36c9SPeter Brune   if (!fas->smoothd) {
1133ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1134ab8d36c9SPeter Brune   }
1135ab8d36c9SPeter Brune   if (!fas->smoothu) {
1136ab8d36c9SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothu);CHKERRQ(ierr);
1137ab8d36c9SPeter Brune   }
1138ab8d36c9SPeter Brune   *smooth = fas->smoothu;
1139ab8d36c9SPeter Brune   PetscFunctionReturn(0);
1140ab8d36c9SPeter Brune }
1141d6ad1212SPeter Brune 
1142d6ad1212SPeter Brune #undef __FUNCT__
1143d6ad1212SPeter Brune #define __FUNCT__ "SNESFASGetCoarseSolve"
1144d6ad1212SPeter Brune /*@
1145d6ad1212SPeter Brune    SNESFASGetCoarseSolve - Gets the coarsest solver.
1146d6ad1212SPeter Brune 
1147d6ad1212SPeter Brune    Input Parameters:
1148d6ad1212SPeter Brune +  snes   - the multigrid context
1149d6ad1212SPeter Brune 
1150d6ad1212SPeter Brune    Output Parameters:
1151d6ad1212SPeter Brune    solve  - the coarse-level solver
1152d6ad1212SPeter Brune 
1153d6ad1212SPeter Brune    Level: advanced
1154d6ad1212SPeter Brune 
1155d6ad1212SPeter Brune .keywords: FAS, MG, get, multigrid, solver, coarse
1156d6ad1212SPeter Brune 
1157d6ad1212SPeter Brune .seealso: SNESFASSetInjection(), SNESFASSetRestriction()
1158d6ad1212SPeter Brune @*/
115922d28d08SBarry Smith PetscErrorCode SNESFASGetCoarseSolve(SNES snes, SNES *smooth)
116022d28d08SBarry Smith {
1161d6ad1212SPeter Brune   SNES_FAS       *fas;
1162d6ad1212SPeter Brune   PetscErrorCode ierr;
1163d6ad1212SPeter Brune   SNES           levelsnes;
116422d28d08SBarry Smith 
1165d6ad1212SPeter Brune   PetscFunctionBegin;
1166d6ad1212SPeter Brune   ierr = SNESFASGetCycleSNES(snes, 0, &levelsnes);CHKERRQ(ierr);
1167d6ad1212SPeter Brune   fas = (SNES_FAS *)levelsnes->data;
1168d6ad1212SPeter Brune   /* if the user chooses to differentiate smoothers, create them both at this point */
1169d6ad1212SPeter Brune   if (!fas->smoothd) {
1170d6ad1212SPeter Brune     ierr = SNESFASCycleCreateSmoother_Private(snes, &fas->smoothd);CHKERRQ(ierr);
1171d6ad1212SPeter Brune   }
1172d6ad1212SPeter Brune   *smooth = fas->smoothd;
1173d6ad1212SPeter Brune   PetscFunctionReturn(0);
1174d6ad1212SPeter Brune }
1175